xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 0fdf79fb08699bf9be0aa4d8ba0185e387a216c8)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/
641b6fd38SMatthew G. Knepley #include <petsc/private/kspimpl.h>
71e25c274SJed Brown #include <petscdm.h>
8391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC, PetscBool *);
94b9ad928SBarry Smith 
10f3b08a26SMatthew G. Knepley /*
11f3b08a26SMatthew G. Knepley    Contains the list of registered coarse space construction routines
12f3b08a26SMatthew G. Knepley */
13f3b08a26SMatthew G. Knepley PetscFunctionList PCMGCoarseList = NULL;
14f3b08a26SMatthew G. Knepley 
15d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGMCycle_Private(PC pc, PC_MG_Levels **mglevelsin, PetscBool transpose, PetscBool matapp, PCRichardsonConvergedReason *reason)
16d71ae5a4SJacob Faibussowitsch {
1731567311SBarry Smith   PC_MG        *mg = (PC_MG *)pc->data;
1831567311SBarry Smith   PC_MG_Levels *mgc, *mglevels = *mglevelsin;
1931567311SBarry Smith   PetscInt      cycles = (mglevels->level == 1) ? 1 : (PetscInt)mglevels->cycles;
204b9ad928SBarry Smith 
214b9ad928SBarry Smith   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve, 0, 0, 0, 0));
23fcb023d4SJed Brown   if (!transpose) {
2430b0564aSStefano Zampini     if (matapp) {
259566063dSJacob Faibussowitsch       PetscCall(KSPMatSolve(mglevels->smoothd, mglevels->B, mglevels->X)); /* pre-smooth */
269566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->smoothd, pc, NULL));
2730b0564aSStefano Zampini     } else {
289566063dSJacob Faibussowitsch       PetscCall(KSPSolve(mglevels->smoothd, mglevels->b, mglevels->x)); /* pre-smooth */
299566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->smoothd, pc, mglevels->x));
3030b0564aSStefano Zampini     }
31fcb023d4SJed Brown   } else {
3228b400f6SJacob Faibussowitsch     PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported");
339566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(mglevels->smoothu, mglevels->b, mglevels->x)); /* transpose of post-smooth */
349566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(mglevels->smoothu, pc, mglevels->x));
35fcb023d4SJed Brown   }
369566063dSJacob Faibussowitsch   if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve, 0, 0, 0, 0));
3731567311SBarry Smith   if (mglevels->level) { /* not the coarsest grid */
389566063dSJacob Faibussowitsch     if (mglevels->eventresidual) PetscCall(PetscLogEventBegin(mglevels->eventresidual, 0, 0, 0, 0));
3948a46eb9SPierre Jolivet     if (matapp && !mglevels->R) PetscCall(MatDuplicate(mglevels->B, MAT_DO_NOT_COPY_VALUES, &mglevels->R));
40fcb023d4SJed Brown     if (!transpose) {
419566063dSJacob Faibussowitsch       if (matapp) PetscCall((*mglevels->matresidual)(mglevels->A, mglevels->B, mglevels->X, mglevels->R));
429566063dSJacob Faibussowitsch       else PetscCall((*mglevels->residual)(mglevels->A, mglevels->b, mglevels->x, mglevels->r));
43fcb023d4SJed Brown     } else {
449566063dSJacob Faibussowitsch       if (matapp) PetscCall((*mglevels->matresidualtranspose)(mglevels->A, mglevels->B, mglevels->X, mglevels->R));
459566063dSJacob Faibussowitsch       else PetscCall((*mglevels->residualtranspose)(mglevels->A, mglevels->b, mglevels->x, mglevels->r));
46fcb023d4SJed Brown     }
479566063dSJacob Faibussowitsch     if (mglevels->eventresidual) PetscCall(PetscLogEventEnd(mglevels->eventresidual, 0, 0, 0, 0));
484b9ad928SBarry Smith 
494b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
5031567311SBarry Smith     if (mglevels->level == mglevels->levels - 1 && mg->ttol && reason) {
514b9ad928SBarry Smith       PetscReal rnorm;
529566063dSJacob Faibussowitsch       PetscCall(VecNorm(mglevels->r, NORM_2, &rnorm));
534b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
5470441072SBarry Smith         if (rnorm < mg->abstol) {
554d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
569566063dSJacob Faibussowitsch           PetscCall(PetscInfo(pc, "Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n", (double)rnorm, (double)mg->abstol));
574b9ad928SBarry Smith         } else {
584d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
599566063dSJacob 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));
604b9ad928SBarry Smith         }
614b9ad928SBarry Smith         PetscFunctionReturn(0);
624b9ad928SBarry Smith       }
634b9ad928SBarry Smith     }
644b9ad928SBarry Smith 
6531567311SBarry Smith     mgc = *(mglevelsin - 1);
669566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict, 0, 0, 0, 0));
67fcb023d4SJed Brown     if (!transpose) {
689566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatRestrict(mglevels->restrct, mglevels->R, &mgc->B));
699566063dSJacob Faibussowitsch       else PetscCall(MatRestrict(mglevels->restrct, mglevels->r, mgc->b));
70fcb023d4SJed Brown     } else {
719566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatRestrict(mglevels->interpolate, mglevels->R, &mgc->B));
729566063dSJacob Faibussowitsch       else PetscCall(MatRestrict(mglevels->interpolate, mglevels->r, mgc->b));
73fcb023d4SJed Brown     }
749566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict, 0, 0, 0, 0));
7530b0564aSStefano Zampini     if (matapp) {
7630b0564aSStefano Zampini       if (!mgc->X) {
779566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(mgc->B, MAT_DO_NOT_COPY_VALUES, &mgc->X));
7830b0564aSStefano Zampini       } else {
799566063dSJacob Faibussowitsch         PetscCall(MatZeroEntries(mgc->X));
8030b0564aSStefano Zampini       }
8130b0564aSStefano Zampini     } else {
829566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(mgc->x));
8330b0564aSStefano Zampini     }
8448a46eb9SPierre Jolivet     while (cycles--) PetscCall(PCMGMCycle_Private(pc, mglevelsin - 1, transpose, matapp, reason));
859566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict, 0, 0, 0, 0));
86fcb023d4SJed Brown     if (!transpose) {
879566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatInterpolateAdd(mglevels->interpolate, mgc->X, mglevels->X, &mglevels->X));
889566063dSJacob Faibussowitsch       else PetscCall(MatInterpolateAdd(mglevels->interpolate, mgc->x, mglevels->x, mglevels->x));
89fcb023d4SJed Brown     } else {
909566063dSJacob Faibussowitsch       PetscCall(MatInterpolateAdd(mglevels->restrct, mgc->x, mglevels->x, mglevels->x));
91fcb023d4SJed Brown     }
929566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict, 0, 0, 0, 0));
939566063dSJacob Faibussowitsch     if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve, 0, 0, 0, 0));
94fcb023d4SJed Brown     if (!transpose) {
9530b0564aSStefano Zampini       if (matapp) {
969566063dSJacob Faibussowitsch         PetscCall(KSPMatSolve(mglevels->smoothu, mglevels->B, mglevels->X)); /* post smooth */
979566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(mglevels->smoothu, pc, NULL));
9830b0564aSStefano Zampini       } else {
999566063dSJacob Faibussowitsch         PetscCall(KSPSolve(mglevels->smoothu, mglevels->b, mglevels->x)); /* post smooth */
1009566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(mglevels->smoothu, pc, mglevels->x));
10130b0564aSStefano Zampini       }
102fcb023d4SJed Brown     } else {
10328b400f6SJacob Faibussowitsch       PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported");
1049566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(mglevels->smoothd, mglevels->b, mglevels->x)); /* post smooth */
1059566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->smoothd, pc, mglevels->x));
106fcb023d4SJed Brown     }
10741b6fd38SMatthew G. Knepley     if (mglevels->cr) {
10828b400f6SJacob Faibussowitsch       PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported");
10941b6fd38SMatthew G. Knepley       /* TODO Turn on copy and turn off noisy if we have an exact solution
1109566063dSJacob Faibussowitsch       PetscCall(VecCopy(mglevels->x, mglevels->crx));
1119566063dSJacob Faibussowitsch       PetscCall(VecCopy(mglevels->b, mglevels->crb)); */
1129566063dSJacob Faibussowitsch       PetscCall(KSPSetNoisy_Private(mglevels->crx));
1139566063dSJacob Faibussowitsch       PetscCall(KSPSolve(mglevels->cr, mglevels->crb, mglevels->crx)); /* compatible relaxation */
1149566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->cr, pc, mglevels->crx));
11541b6fd38SMatthew G. Knepley     }
1169566063dSJacob Faibussowitsch     if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve, 0, 0, 0, 0));
1174b9ad928SBarry Smith   }
1184b9ad928SBarry Smith   PetscFunctionReturn(0);
1194b9ad928SBarry Smith }
1204b9ad928SBarry Smith 
121d71ae5a4SJacob 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)
122d71ae5a4SJacob Faibussowitsch {
123f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
124f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
125391689abSStefano Zampini   PC             tpc;
126391689abSStefano Zampini   PetscBool      changeu, changed;
127f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels, i;
1284b9ad928SBarry Smith 
1294b9ad928SBarry Smith   PetscFunctionBegin;
130a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
131a762d673SBarry Smith   for (i = 0; i < levels; i++) {
132a762d673SBarry Smith     if (!mglevels[i]->A) {
1339566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL));
1349566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A));
135a762d673SBarry Smith     }
136a762d673SBarry Smith   }
137391689abSStefano Zampini 
1389566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc));
1399566063dSJacob Faibussowitsch   PetscCall(PCPreSolveChangeRHS(tpc, &changed));
1409566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc));
1419566063dSJacob Faibussowitsch   PetscCall(PCPreSolveChangeRHS(tpc, &changeu));
142391689abSStefano Zampini   if (!changed && !changeu) {
1439566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[levels - 1]->b));
144f3fbd535SBarry Smith     mglevels[levels - 1]->b = b;
145391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
146391689abSStefano Zampini     if (!mglevels[levels - 1]->b) {
147391689abSStefano Zampini       Vec *vec;
148391689abSStefano Zampini 
1499566063dSJacob Faibussowitsch       PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL));
150391689abSStefano Zampini       mglevels[levels - 1]->b = *vec;
1519566063dSJacob Faibussowitsch       PetscCall(PetscFree(vec));
152391689abSStefano Zampini     }
1539566063dSJacob Faibussowitsch     PetscCall(VecCopy(b, mglevels[levels - 1]->b));
154391689abSStefano Zampini   }
155f3fbd535SBarry Smith   mglevels[levels - 1]->x = x;
1564b9ad928SBarry Smith 
15731567311SBarry Smith   mg->rtol   = rtol;
15831567311SBarry Smith   mg->abstol = abstol;
15931567311SBarry Smith   mg->dtol   = dtol;
1604b9ad928SBarry Smith   if (rtol) {
1614b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
1624b9ad928SBarry Smith     PetscReal rnorm;
1637319c654SBarry Smith     if (zeroguess) {
1649566063dSJacob Faibussowitsch       PetscCall(VecNorm(b, NORM_2, &rnorm));
1657319c654SBarry Smith     } else {
1669566063dSJacob Faibussowitsch       PetscCall((*mglevels[levels - 1]->residual)(mglevels[levels - 1]->A, b, x, w));
1679566063dSJacob Faibussowitsch       PetscCall(VecNorm(w, NORM_2, &rnorm));
1687319c654SBarry Smith     }
16931567311SBarry Smith     mg->ttol = PetscMax(rtol * rnorm, abstol);
1702fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
1712fa5cd67SKarl Rupp   else mg->ttol = 0.0;
1724b9ad928SBarry Smith 
1734d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
174336babb1SJed Brown      stop prematurely due to small residual */
1754d0a8057SBarry Smith   for (i = 1; i < levels; i++) {
1769566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[i]->smoothu, 0, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
177f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
17823067569SBarry Smith       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
1799566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE));
1809566063dSJacob Faibussowitsch       PetscCall(KSPSetTolerances(mglevels[i]->smoothd, 0, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
1814b9ad928SBarry Smith     }
1824d0a8057SBarry Smith   }
1834d0a8057SBarry Smith 
1844d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1854d0a8057SBarry Smith   for (i = 0; i < its; i++) {
1869566063dSJacob Faibussowitsch     PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, PETSC_FALSE, PETSC_FALSE, reason));
1874d0a8057SBarry Smith     if (*reason) break;
1884d0a8057SBarry Smith   }
1894d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1904d0a8057SBarry Smith   *outits = i;
191391689abSStefano Zampini   if (!changed && !changeu) mglevels[levels - 1]->b = NULL;
1924b9ad928SBarry Smith   PetscFunctionReturn(0);
1934b9ad928SBarry Smith }
1944b9ad928SBarry Smith 
195d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset_MG(PC pc)
196d71ae5a4SJacob Faibussowitsch {
1973aeaf226SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
1983aeaf226SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
1992b3cbbdaSStefano Zampini   PetscInt       i, n;
2003aeaf226SBarry Smith 
2013aeaf226SBarry Smith   PetscFunctionBegin;
2023aeaf226SBarry Smith   if (mglevels) {
2033aeaf226SBarry Smith     n = mglevels[0]->levels;
2043aeaf226SBarry Smith     for (i = 0; i < n - 1; i++) {
2059566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i + 1]->r));
2069566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->b));
2079566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->x));
2089566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i + 1]->R));
2099566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i]->B));
2109566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i]->X));
2119566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->crx));
2129566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->crb));
2139566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i + 1]->restrct));
2149566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i + 1]->interpolate));
2159566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i + 1]->inject));
2169566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i + 1]->rscale));
2173aeaf226SBarry Smith     }
2189566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[n - 1]->crx));
2199566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[n - 1]->crb));
220391689abSStefano Zampini     /* this is not null only if the smoother on the finest level
221391689abSStefano Zampini        changes the rhs during PreSolve */
2229566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[n - 1]->b));
2239566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&mglevels[n - 1]->B));
2243aeaf226SBarry Smith 
2253aeaf226SBarry Smith     for (i = 0; i < n; i++) {
2262b3cbbdaSStefano Zampini       PetscCall(MatDestroy(&mglevels[i]->coarseSpace));
2279566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i]->A));
22848a46eb9SPierre Jolivet       if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPReset(mglevels[i]->smoothd));
2299566063dSJacob Faibussowitsch       PetscCall(KSPReset(mglevels[i]->smoothu));
2309566063dSJacob Faibussowitsch       if (mglevels[i]->cr) PetscCall(KSPReset(mglevels[i]->cr));
2313aeaf226SBarry Smith     }
232f3b08a26SMatthew G. Knepley     mg->Nc = 0;
2333aeaf226SBarry Smith   }
2343aeaf226SBarry Smith   PetscFunctionReturn(0);
2353aeaf226SBarry Smith }
2363aeaf226SBarry Smith 
23741b6fd38SMatthew G. Knepley /* Implementing CR
23841b6fd38SMatthew G. Knepley 
23941b6fd38SMatthew 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
24041b6fd38SMatthew G. Knepley 
24141b6fd38SMatthew G. Knepley   Inj^T (Inj Inj^T)^{-1} Inj
24241b6fd38SMatthew G. Knepley 
24341b6fd38SMatthew G. Knepley and if Inj is a VecScatter, as it is now in PETSc, we have
24441b6fd38SMatthew G. Knepley 
24541b6fd38SMatthew G. Knepley   Inj^T Inj
24641b6fd38SMatthew G. Knepley 
24741b6fd38SMatthew G. Knepley and
24841b6fd38SMatthew G. Knepley 
24941b6fd38SMatthew G. Knepley   S = I - Inj^T Inj
25041b6fd38SMatthew G. Knepley 
25141b6fd38SMatthew G. Knepley since
25241b6fd38SMatthew G. Knepley 
25341b6fd38SMatthew G. Knepley   Inj S = Inj - (Inj Inj^T) Inj = 0.
25441b6fd38SMatthew G. Knepley 
25541b6fd38SMatthew G. Knepley Brannick suggests
25641b6fd38SMatthew G. Knepley 
25741b6fd38SMatthew G. Knepley   A \to S^T A S  \qquad\mathrm{and}\qquad M \to S^T M S
25841b6fd38SMatthew G. Knepley 
25941b6fd38SMatthew 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
26041b6fd38SMatthew G. Knepley 
26141b6fd38SMatthew G. Knepley   M^{-1} A \to S M^{-1} A S
26241b6fd38SMatthew G. Knepley 
26341b6fd38SMatthew G. Knepley In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left.
26441b6fd38SMatthew G. Knepley 
26541b6fd38SMatthew G. Knepley   Check: || Inj P - I ||_F < tol
26641b6fd38SMatthew G. Knepley   Check: In general, Inj Inj^T = I
26741b6fd38SMatthew G. Knepley */
26841b6fd38SMatthew G. Knepley 
26941b6fd38SMatthew G. Knepley typedef struct {
27041b6fd38SMatthew G. Knepley   PC       mg;  /* The PCMG object */
27141b6fd38SMatthew G. Knepley   PetscInt l;   /* The multigrid level for this solver */
27241b6fd38SMatthew G. Knepley   Mat      Inj; /* The injection matrix */
27341b6fd38SMatthew G. Knepley   Mat      S;   /* I - Inj^T Inj */
27441b6fd38SMatthew G. Knepley } CRContext;
27541b6fd38SMatthew G. Knepley 
276d71ae5a4SJacob Faibussowitsch static PetscErrorCode CRSetup_Private(PC pc)
277d71ae5a4SJacob Faibussowitsch {
27841b6fd38SMatthew G. Knepley   CRContext *ctx;
27941b6fd38SMatthew G. Knepley   Mat        It;
28041b6fd38SMatthew G. Knepley 
28141b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
2829566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &ctx));
2839566063dSJacob Faibussowitsch   PetscCall(PCMGGetInjection(ctx->mg, ctx->l, &It));
28428b400f6SJacob Faibussowitsch   PetscCheck(It, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG");
2859566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(It, &ctx->Inj));
2869566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(ctx->Inj, &ctx->S));
2879566063dSJacob Faibussowitsch   PetscCall(MatScale(ctx->S, -1.0));
2889566063dSJacob Faibussowitsch   PetscCall(MatShift(ctx->S, 1.0));
28941b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
29041b6fd38SMatthew G. Knepley }
29141b6fd38SMatthew G. Knepley 
292d71ae5a4SJacob Faibussowitsch static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y)
293d71ae5a4SJacob Faibussowitsch {
29441b6fd38SMatthew G. Knepley   CRContext *ctx;
29541b6fd38SMatthew G. Knepley 
29641b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
2979566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &ctx));
2989566063dSJacob Faibussowitsch   PetscCall(MatMult(ctx->S, x, y));
29941b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
30041b6fd38SMatthew G. Knepley }
30141b6fd38SMatthew G. Knepley 
302d71ae5a4SJacob Faibussowitsch static PetscErrorCode CRDestroy_Private(PC pc)
303d71ae5a4SJacob Faibussowitsch {
30441b6fd38SMatthew G. Knepley   CRContext *ctx;
30541b6fd38SMatthew G. Knepley 
30641b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
3079566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &ctx));
3089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->Inj));
3099566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->S));
3109566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
3119566063dSJacob Faibussowitsch   PetscCall(PCShellSetContext(pc, NULL));
31241b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
31341b6fd38SMatthew G. Knepley }
31441b6fd38SMatthew G. Knepley 
315d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr)
316d71ae5a4SJacob Faibussowitsch {
31741b6fd38SMatthew G. Knepley   CRContext *ctx;
31841b6fd38SMatthew G. Knepley 
31941b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
3209566063dSJacob Faibussowitsch   PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), cr));
3219566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*cr, "S (complementary projector to injection)"));
3229566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(1, &ctx));
32341b6fd38SMatthew G. Knepley   ctx->mg = pc;
32441b6fd38SMatthew G. Knepley   ctx->l  = l;
3259566063dSJacob Faibussowitsch   PetscCall(PCSetType(*cr, PCSHELL));
3269566063dSJacob Faibussowitsch   PetscCall(PCShellSetContext(*cr, ctx));
3279566063dSJacob Faibussowitsch   PetscCall(PCShellSetApply(*cr, CRApply_Private));
3289566063dSJacob Faibussowitsch   PetscCall(PCShellSetSetUp(*cr, CRSetup_Private));
3299566063dSJacob Faibussowitsch   PetscCall(PCShellSetDestroy(*cr, CRDestroy_Private));
33041b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
33141b6fd38SMatthew G. Knepley }
33241b6fd38SMatthew G. Knepley 
333d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetLevels_MG(PC pc, PetscInt levels, MPI_Comm *comms)
334d71ae5a4SJacob Faibussowitsch {
335f3fbd535SBarry Smith   PC_MG         *mg = (PC_MG *)pc->data;
336ce94432eSBarry Smith   MPI_Comm       comm;
3373aeaf226SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
33810eca3edSLisandro Dalcin   PCMGType       mgtype   = mg->am;
33910167fecSBarry Smith   PetscInt       mgctype  = (PetscInt)PC_MG_CYCLE_V;
340f3fbd535SBarry Smith   PetscInt       i;
341f3fbd535SBarry Smith   PetscMPIInt    size;
342f3fbd535SBarry Smith   const char    *prefix;
343f3fbd535SBarry Smith   PC             ipc;
3443aeaf226SBarry Smith   PetscInt       n;
3454b9ad928SBarry Smith 
3464b9ad928SBarry Smith   PetscFunctionBegin;
3470700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
348548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc, levels, 2);
349548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
3509566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
3513aeaf226SBarry Smith   if (mglevels) {
35210eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
3533aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
3549566063dSJacob Faibussowitsch     PetscCall(PCReset_MG(pc));
3553aeaf226SBarry Smith     n = mglevels[0]->levels;
3563aeaf226SBarry Smith     for (i = 0; i < n; i++) {
35748a46eb9SPierre Jolivet       if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPDestroy(&mglevels[i]->smoothd));
3589566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&mglevels[i]->smoothu));
3599566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&mglevels[i]->cr));
3609566063dSJacob Faibussowitsch       PetscCall(PetscFree(mglevels[i]));
3613aeaf226SBarry Smith     }
3629566063dSJacob Faibussowitsch     PetscCall(PetscFree(mg->levels));
3633aeaf226SBarry Smith   }
364f3fbd535SBarry Smith 
365f3fbd535SBarry Smith   mg->nlevels = levels;
366f3fbd535SBarry Smith 
3679566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(levels, &mglevels));
368f3fbd535SBarry Smith 
3699566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &prefix));
370f3fbd535SBarry Smith 
371a9db87a2SMatthew G Knepley   mg->stageApply = 0;
372f3fbd535SBarry Smith   for (i = 0; i < levels; i++) {
3734dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew(&mglevels[i]));
3742fa5cd67SKarl Rupp 
375f3fbd535SBarry Smith     mglevels[i]->level               = i;
376f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
37710eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
378336babb1SJed Brown     mg->default_smoothu              = 2;
379336babb1SJed Brown     mg->default_smoothd              = 2;
38063e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
38163e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
38263e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
38363e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
384f3fbd535SBarry Smith 
385f3fbd535SBarry Smith     if (comms) comm = comms[i];
386d5a8f7e6SBarry Smith     if (comm != MPI_COMM_NULL) {
3879566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm, &mglevels[i]->smoothd));
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));
4284b9ad928SBarry Smith   PetscFunctionReturn(0);
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 
435f1580f4eSBarry Smith    Logically Collective on pc
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
442f1580f4eSBarry 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:
44897d33e41SMatthew G. Knepley      If the number of levels is one then the multigrid uses the -mg_levels prefix
44997d33e41SMatthew G. Knepley      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
466d5a8f7e6SBarry Smith      recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries.
467d5a8f7e6SBarry Smith 
468f1580f4eSBarry Smith    Fortran Note:
469f1580f4eSBarry 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 
472db781477SPatrick Sanan .seealso: `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);
47897d33e41SMatthew G. Knepley   if (comms) PetscValidPointer(comms, 3);
479cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGSetLevels_C", (PC, PetscInt, MPI_Comm *), (pc, levels, comms));
48097d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
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));
515f3fbd535SBarry Smith   PetscFunctionReturn(0);
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   }
627f3fbd535SBarry Smith   PetscFunctionReturn(0);
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));
634fcb023d4SJed Brown   PetscFunctionReturn(0);
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));
64130b0564aSStefano Zampini   PetscFunctionReturn(0);
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));
648fcb023d4SJed Brown   PetscFunctionReturn(0);
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 *)&gtype, &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 
718ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
719eec35531SMatthew G Knepley     {
720eec35531SMatthew G Knepley       const char   *sname = "MG Apply";
721eec35531SMatthew G Knepley       PetscStageLog stageLog;
722eec35531SMatthew G Knepley       PetscInt      st;
723eec35531SMatthew G Knepley 
7249566063dSJacob Faibussowitsch       PetscCall(PetscLogGetStageLog(&stageLog));
725eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
726eec35531SMatthew G Knepley         PetscBool same;
727eec35531SMatthew G Knepley 
7289566063dSJacob Faibussowitsch         PetscCall(PetscStrcmp(stageLog->stageInfo[st].name, sname, &same));
7292fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
730eec35531SMatthew G Knepley       }
73148a46eb9SPierre Jolivet       if (!mg->stageApply) PetscCall(PetscLogStageRegister(sname, &mg->stageApply));
732eec35531SMatthew G Knepley     }
733ec60ca73SSatish Balay #endif
734f3fbd535SBarry Smith   }
735d0609cedSBarry Smith   PetscOptionsHeadEnd();
736f3b08a26SMatthew G. Knepley   /* Check option consistency */
7379566063dSJacob Faibussowitsch   PetscCall(PCMGGetGalerkin(pc, &gtype));
7389566063dSJacob Faibussowitsch   PetscCall(PCMGGetAdaptInterpolation(pc, &flg));
73908401ef6SPierre Jolivet   PetscCheck(!flg || !(gtype >= PC_MG_GALERKIN_NONE), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
740f3fbd535SBarry Smith   PetscFunctionReturn(0);
741f3fbd535SBarry Smith }
742f3fbd535SBarry Smith 
7430a545947SLisandro Dalcin const char *const PCMGTypes[]            = {"MULTIPLICATIVE", "ADDITIVE", "FULL", "KASKADE", "PCMGType", "PC_MG", NULL};
7440a545947SLisandro Dalcin const char *const PCMGCycleTypes[]       = {"invalid", "v", "w", "PCMGCycleType", "PC_MG_CYCLE", NULL};
7450a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[]    = {"both", "pmat", "mat", "none", "external", "PCMGGalerkinType", "PC_MG_GALERKIN", NULL};
7462b3cbbdaSStefano Zampini const char *const PCMGCoarseSpaceTypes[] = {"none", "polynomial", "harmonic", "eigenvector", "generalized_eigenvector", "gdsw", "PCMGCoarseSpaceType", "PCMG_ADAPT_NONE", NULL};
747f3fbd535SBarry Smith 
7489804daf3SBarry Smith #include <petscdraw.h>
749d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView_MG(PC pc, PetscViewer viewer)
750d71ae5a4SJacob Faibussowitsch {
751f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
752f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
753e3deeeafSJed Brown   PetscInt       levels   = mglevels ? mglevels[0]->levels : 0, i;
754219b1045SBarry Smith   PetscBool      iascii, isbinary, isdraw;
755f3fbd535SBarry Smith 
756f3fbd535SBarry Smith   PetscFunctionBegin;
7579566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
7589566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
7599566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
760f3fbd535SBarry Smith   if (iascii) {
761e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
76263a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type is %s, levels=%" PetscInt_FMT " cycles=%s\n", PCMGTypes[mg->am], levels, cyclename));
76348a46eb9SPierre Jolivet     if (mg->am == PC_MG_MULTIPLICATIVE) PetscCall(PetscViewerASCIIPrintf(viewer, "    Cycles per PCApply=%" PetscInt_FMT "\n", mg->cyclesperpcapply));
7642134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
7659566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using Galerkin computed coarse grid matrices\n"));
7662134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
7679566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using Galerkin computed coarse grid matrices for pmat\n"));
7682134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
7699566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using Galerkin computed coarse grid matrices for mat\n"));
7702134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
7719566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using externally compute Galerkin coarse grid matrices\n"));
7724f66f45eSBarry Smith     } else {
7739566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Not using Galerkin computed coarse grid matrices\n"));
774f3fbd535SBarry Smith     }
7751baa6e33SBarry Smith     if (mg->view) PetscCall((*mg->view)(pc, viewer));
776f3fbd535SBarry Smith     for (i = 0; i < levels; i++) {
77763a3b9bcSJacob Faibussowitsch       if (i) {
77863a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n", i));
779f3fbd535SBarry Smith       } else {
78063a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n", i));
781f3fbd535SBarry Smith       }
7829566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
7839566063dSJacob Faibussowitsch       PetscCall(KSPView(mglevels[i]->smoothd, viewer));
7849566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
785f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
7869566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) same as down solver (pre-smoother)\n"));
787f3fbd535SBarry Smith       } else if (i) {
78863a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n", i));
7899566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
7909566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothu, viewer));
7919566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
792f3fbd535SBarry Smith       }
79341b6fd38SMatthew G. Knepley       if (i && mglevels[i]->cr) {
79463a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "CR solver on level %" PetscInt_FMT " -------------------------------\n", i));
7959566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
7969566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->cr, viewer));
7979566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
79841b6fd38SMatthew G. Knepley       }
799f3fbd535SBarry Smith     }
8005b0b0462SBarry Smith   } else if (isbinary) {
8015b0b0462SBarry Smith     for (i = levels - 1; i >= 0; i--) {
8029566063dSJacob Faibussowitsch       PetscCall(KSPView(mglevels[i]->smoothd, viewer));
80348a46eb9SPierre Jolivet       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPView(mglevels[i]->smoothu, viewer));
8045b0b0462SBarry Smith     }
805219b1045SBarry Smith   } else if (isdraw) {
806219b1045SBarry Smith     PetscDraw draw;
807b4375e8dSPeter Brune     PetscReal x, w, y, bottom, th;
8089566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
8099566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
8109566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringGetSize(draw, NULL, &th));
811219b1045SBarry Smith     bottom = y - th;
812219b1045SBarry Smith     for (i = levels - 1; i >= 0; i--) {
813b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
8149566063dSJacob Faibussowitsch         PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
8159566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothd, viewer));
8169566063dSJacob Faibussowitsch         PetscCall(PetscDrawPopCurrentPoint(draw));
817b4375e8dSPeter Brune       } else {
818b4375e8dSPeter Brune         w = 0.5 * PetscMin(1.0 - x, x);
8199566063dSJacob Faibussowitsch         PetscCall(PetscDrawPushCurrentPoint(draw, x + w, bottom));
8209566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothd, viewer));
8219566063dSJacob Faibussowitsch         PetscCall(PetscDrawPopCurrentPoint(draw));
8229566063dSJacob Faibussowitsch         PetscCall(PetscDrawPushCurrentPoint(draw, x - w, bottom));
8239566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothu, viewer));
8249566063dSJacob Faibussowitsch         PetscCall(PetscDrawPopCurrentPoint(draw));
825b4375e8dSPeter Brune       }
8269566063dSJacob Faibussowitsch       PetscCall(PetscDrawGetBoundingBox(draw, NULL, &bottom, NULL, NULL));
8271cd381d6SBarry Smith       bottom -= th;
828219b1045SBarry Smith     }
8295b0b0462SBarry Smith   }
830f3fbd535SBarry Smith   PetscFunctionReturn(0);
831f3fbd535SBarry Smith }
832f3fbd535SBarry Smith 
833af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
834cab2e9ccSBarry Smith 
835f3fbd535SBarry Smith /*
836f3fbd535SBarry Smith     Calls setup for the KSP on each level
837f3fbd535SBarry Smith */
838d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp_MG(PC pc)
839d71ae5a4SJacob Faibussowitsch {
840f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
841f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
8421a97d7b6SStefano Zampini   PetscInt       i, n;
84398e04984SBarry Smith   PC             cpc;
8443db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE, opsset, use_amat, missinginterpolate = PETSC_FALSE;
845f3fbd535SBarry Smith   Mat            dA, dB;
846f3fbd535SBarry Smith   Vec            tvec;
847218a07d4SBarry Smith   DM            *dms;
8480a545947SLisandro Dalcin   PetscViewer    viewer = NULL;
84941b6fd38SMatthew G. Knepley   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE;
8502b3cbbdaSStefano Zampini   PetscBool      adaptInterpolation = mg->adaptInterpolation;
851f3fbd535SBarry Smith 
852f3fbd535SBarry Smith   PetscFunctionBegin;
85328b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels with PCMGSetLevels() before setting up");
8541a97d7b6SStefano Zampini   n = mglevels[0]->levels;
85501bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
8563aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
8573aeaf226SBarry Smith     PetscInt levels;
8589566063dSJacob Faibussowitsch     PetscCall(DMGetRefineLevel(pc->dm, &levels));
8593aeaf226SBarry Smith     levels++;
8603aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
8619566063dSJacob Faibussowitsch       PetscCall(PCMGSetLevels(pc, levels, NULL));
8623aeaf226SBarry Smith       n = levels;
8639566063dSJacob Faibussowitsch       PetscCall(PCSetFromOptions(pc)); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
8643aeaf226SBarry Smith       mglevels = mg->levels;
8653aeaf226SBarry Smith     }
8663aeaf226SBarry Smith   }
8679566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[0]->smoothd, &cpc));
8683aeaf226SBarry Smith 
869f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
870f3fbd535SBarry Smith   /* so use those from global PC */
871f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
8729566063dSJacob Faibussowitsch   PetscCall(KSPGetOperatorsSet(mglevels[n - 1]->smoothd, NULL, &opsset));
87398e04984SBarry Smith   if (opsset) {
87498e04984SBarry Smith     Mat mmat;
8759566063dSJacob Faibussowitsch     PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, NULL, &mmat));
87698e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
87798e04984SBarry Smith   }
8784a5f07a5SMark F. Adams 
87941b6fd38SMatthew G. Knepley   /* Create CR solvers */
8809566063dSJacob Faibussowitsch   PetscCall(PCMGGetAdaptCR(pc, &doCR));
88141b6fd38SMatthew G. Knepley   if (doCR) {
88241b6fd38SMatthew G. Knepley     const char *prefix;
88341b6fd38SMatthew G. Knepley 
8849566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
88541b6fd38SMatthew G. Knepley     for (i = 1; i < n; ++i) {
88641b6fd38SMatthew G. Knepley       PC   ipc, cr;
88741b6fd38SMatthew G. Knepley       char crprefix[128];
88841b6fd38SMatthew G. Knepley 
8899566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &mglevels[i]->cr));
8909566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE));
8919566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->cr, (PetscObject)pc, n - i));
8929566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(mglevels[i]->cr, prefix));
8939566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level));
8949566063dSJacob Faibussowitsch       PetscCall(KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV));
8959566063dSJacob Faibussowitsch       PetscCall(KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL));
8969566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED));
8979566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(mglevels[i]->cr, &ipc));
89841b6fd38SMatthew G. Knepley 
8999566063dSJacob Faibussowitsch       PetscCall(PCSetType(ipc, PCCOMPOSITE));
9009566063dSJacob Faibussowitsch       PetscCall(PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE));
9019566063dSJacob Faibussowitsch       PetscCall(PCCompositeAddPCType(ipc, PCSOR));
9029566063dSJacob Faibussowitsch       PetscCall(CreateCR_Private(pc, i, &cr));
9039566063dSJacob Faibussowitsch       PetscCall(PCCompositeAddPC(ipc, cr));
9049566063dSJacob Faibussowitsch       PetscCall(PCDestroy(&cr));
90541b6fd38SMatthew G. Knepley 
9069566063dSJacob Faibussowitsch       PetscCall(KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd));
9079566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE));
9089566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int)i));
9099566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix));
91041b6fd38SMatthew G. Knepley     }
91141b6fd38SMatthew G. Knepley   }
91241b6fd38SMatthew G. Knepley 
9134a5f07a5SMark F. Adams   if (!opsset) {
9149566063dSJacob Faibussowitsch     PetscCall(PCGetUseAmat(pc, &use_amat));
9154a5f07a5SMark F. Adams     if (use_amat) {
9169566063dSJacob 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"));
9179566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->mat, pc->pmat));
91869858f1bSStefano Zampini     } else {
9199566063dSJacob 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"));
9209566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->pmat, pc->pmat));
9214a5f07a5SMark F. Adams     }
9224a5f07a5SMark F. Adams   }
923f3fbd535SBarry Smith 
9249c7ffb39SBarry Smith   for (i = n - 1; i > 0; i--) {
9259c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
9269c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
9272b3cbbdaSStefano Zampini       break;
9289c7ffb39SBarry Smith     }
9299c7ffb39SBarry Smith   }
9302134b1e4SBarry Smith 
9319566063dSJacob Faibussowitsch   PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, &dA, &dB));
9322134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
9332b3cbbdaSStefano Zampini   if (mg->galerkin == PC_MG_GALERKIN_NONE || ((mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_MAT) && !dAeqdB)) {
9342134b1e4SBarry Smith     needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
9352134b1e4SBarry Smith   }
9362134b1e4SBarry Smith 
9372b3cbbdaSStefano Zampini   if (pc->dm && !pc->setupcalled) {
9382b3cbbdaSStefano Zampini     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
9392b3cbbdaSStefano Zampini     PetscCall(KSPSetDM(mglevels[n - 1]->smoothd, pc->dm));
9402b3cbbdaSStefano Zampini     PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothd, PETSC_FALSE));
9412b3cbbdaSStefano Zampini     if (mglevels[n - 1]->smoothd != mglevels[n - 1]->smoothu) {
9422b3cbbdaSStefano Zampini       PetscCall(KSPSetDM(mglevels[n - 1]->smoothu, pc->dm));
9432b3cbbdaSStefano Zampini       PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothu, PETSC_FALSE));
9442b3cbbdaSStefano Zampini     }
9452b3cbbdaSStefano Zampini     if (mglevels[n - 1]->cr) {
9462b3cbbdaSStefano Zampini       PetscCall(KSPSetDM(mglevels[n - 1]->cr, pc->dm));
9472b3cbbdaSStefano Zampini       PetscCall(KSPSetDMActive(mglevels[n - 1]->cr, PETSC_FALSE));
9482b3cbbdaSStefano Zampini     }
9492b3cbbdaSStefano Zampini   }
9502b3cbbdaSStefano Zampini 
9519c7ffb39SBarry Smith   /*
9529c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
9532b3cbbdaSStefano Zampini    Skipping for externally managed hierarchy (such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
9549c7ffb39SBarry Smith   */
95532fe681dSStefano Zampini   if (missinginterpolate && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
95632fe681dSStefano Zampini     /* first see if we can compute a coarse space */
95732fe681dSStefano Zampini     if (mg->coarseSpaceType == PCMG_ADAPT_GDSW) {
95832fe681dSStefano Zampini       for (i = n - 2; i > -1; i--) {
95932fe681dSStefano Zampini         if (!mglevels[i + 1]->restrct && !mglevels[i + 1]->interpolate) {
96032fe681dSStefano Zampini           PetscCall(PCMGComputeCoarseSpace_Internal(pc, i + 1, mg->coarseSpaceType, mg->Nc, NULL, &mglevels[i + 1]->coarseSpace));
96132fe681dSStefano Zampini           PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->coarseSpace));
96232fe681dSStefano Zampini         }
96332fe681dSStefano Zampini       }
96432fe681dSStefano Zampini     } else { /* construct the interpolation from the DMs */
9652e499ae9SBarry Smith       Mat p;
96673250ac0SBarry Smith       Vec rscale;
9679566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &dms));
968218a07d4SBarry Smith       dms[n - 1] = pc->dm;
969ef1267afSMatthew G. Knepley       /* Separately create them so we do not get DMKSP interference between levels */
9709566063dSJacob Faibussowitsch       for (i = n - 2; i > -1; i--) PetscCall(DMCoarsen(dms[i + 1], MPI_COMM_NULL, &dms[i]));
971218a07d4SBarry Smith       for (i = n - 2; i > -1; i--) {
972942e3340SBarry Smith         DMKSP     kdm;
973eab52d2bSLawrence Mitchell         PetscBool dmhasrestrict, dmhasinject;
9749566063dSJacob Faibussowitsch         PetscCall(KSPSetDM(mglevels[i]->smoothd, dms[i]));
9759566063dSJacob Faibussowitsch         if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothd, PETSC_FALSE));
976c27ee7a3SPatrick Farrell         if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
9779566063dSJacob Faibussowitsch           PetscCall(KSPSetDM(mglevels[i]->smoothu, dms[i]));
9789566063dSJacob Faibussowitsch           if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothu, PETSC_FALSE));
979c27ee7a3SPatrick Farrell         }
98041b6fd38SMatthew G. Knepley         if (mglevels[i]->cr) {
9819566063dSJacob Faibussowitsch           PetscCall(KSPSetDM(mglevels[i]->cr, dms[i]));
9829566063dSJacob Faibussowitsch           if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->cr, PETSC_FALSE));
98341b6fd38SMatthew G. Knepley         }
9849566063dSJacob Faibussowitsch         PetscCall(DMGetDMKSPWrite(dms[i], &kdm));
985d1a3e677SJed Brown         /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
986d1a3e677SJed Brown          * a bitwise OR of computing the matrix, RHS, and initial iterate. */
9870298fd71SBarry Smith         kdm->ops->computerhs = NULL;
9880298fd71SBarry Smith         kdm->rhsctx          = NULL;
98924c3aa18SBarry Smith         if (!mglevels[i + 1]->interpolate) {
9909566063dSJacob Faibussowitsch           PetscCall(DMCreateInterpolation(dms[i], dms[i + 1], &p, &rscale));
9919566063dSJacob Faibussowitsch           PetscCall(PCMGSetInterpolation(pc, i + 1, p));
9929566063dSJacob Faibussowitsch           if (rscale) PetscCall(PCMGSetRScale(pc, i + 1, rscale));
9939566063dSJacob Faibussowitsch           PetscCall(VecDestroy(&rscale));
9949566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&p));
995218a07d4SBarry Smith         }
9969566063dSJacob Faibussowitsch         PetscCall(DMHasCreateRestriction(dms[i], &dmhasrestrict));
9973ad4599aSBarry Smith         if (dmhasrestrict && !mglevels[i + 1]->restrct) {
9989566063dSJacob Faibussowitsch           PetscCall(DMCreateRestriction(dms[i], dms[i + 1], &p));
9999566063dSJacob Faibussowitsch           PetscCall(PCMGSetRestriction(pc, i + 1, p));
10009566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&p));
10013ad4599aSBarry Smith         }
10029566063dSJacob Faibussowitsch         PetscCall(DMHasCreateInjection(dms[i], &dmhasinject));
1003eab52d2bSLawrence Mitchell         if (dmhasinject && !mglevels[i + 1]->inject) {
10049566063dSJacob Faibussowitsch           PetscCall(DMCreateInjection(dms[i], dms[i + 1], &p));
10059566063dSJacob Faibussowitsch           PetscCall(PCMGSetInjection(pc, i + 1, p));
10069566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&p));
1007eab52d2bSLawrence Mitchell         }
100824c3aa18SBarry Smith       }
10092d2b81a6SBarry Smith 
10109566063dSJacob Faibussowitsch       for (i = n - 2; i > -1; i--) PetscCall(DMDestroy(&dms[i]));
10119566063dSJacob Faibussowitsch       PetscCall(PetscFree(dms));
1012ce4cda84SJed Brown     }
101332fe681dSStefano Zampini   }
1014cab2e9ccSBarry Smith 
10152134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
10162134b1e4SBarry Smith     Mat       A, B;
10172134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE, doB = PETSC_FALSE;
10182134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
10192134b1e4SBarry Smith 
10202b3cbbdaSStefano Zampini     if (mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_BOTH) doB = PETSC_TRUE;
10212b3cbbdaSStefano Zampini     if (mg->galerkin == PC_MG_GALERKIN_MAT || (mg->galerkin == PC_MG_GALERKIN_BOTH && dA != dB)) doA = PETSC_TRUE;
10222134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
1023f3fbd535SBarry Smith     for (i = n - 2; i > -1; i--) {
10247827d75bSBarry 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");
102548a46eb9SPierre Jolivet       if (!mglevels[i + 1]->interpolate) PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->restrct));
102648a46eb9SPierre Jolivet       if (!mglevels[i + 1]->restrct) PetscCall(PCMGSetRestriction(pc, i + 1, mglevels[i + 1]->interpolate));
102748a46eb9SPierre Jolivet       if (reuse == MAT_REUSE_MATRIX) PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, &B));
102848a46eb9SPierre Jolivet       if (doA) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dA, mglevels[i + 1]->interpolate, reuse, 1.0, &A));
102948a46eb9SPierre Jolivet       if (doB) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dB, mglevels[i + 1]->interpolate, reuse, 1.0, &B));
10302134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
10312134b1e4SBarry Smith       if (!doA && dAeqdB) {
10329566063dSJacob Faibussowitsch         if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)B));
10332134b1e4SBarry Smith         A = B;
10342134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
10359566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, NULL));
10369566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)A));
1037b9367914SBarry Smith       }
10382134b1e4SBarry Smith       if (!doB && dAeqdB) {
10399566063dSJacob Faibussowitsch         if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)A));
10402134b1e4SBarry Smith         B = A;
10412134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
10429566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd, NULL, &B));
10439566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)B));
10447d537d92SJed Brown       }
10452134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
10469566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[i]->smoothd, A, B));
10479566063dSJacob Faibussowitsch         PetscCall(PetscObjectDereference((PetscObject)A));
10489566063dSJacob Faibussowitsch         PetscCall(PetscObjectDereference((PetscObject)B));
10492134b1e4SBarry Smith       }
10502134b1e4SBarry Smith       dA = A;
1051f3fbd535SBarry Smith       dB = B;
1052f3fbd535SBarry Smith     }
1053f3fbd535SBarry Smith   }
1054f3b08a26SMatthew G. Knepley 
1055f3b08a26SMatthew G. Knepley   /* Adapt interpolation matrices */
10562b3cbbdaSStefano Zampini   if (adaptInterpolation) {
1057f3b08a26SMatthew G. Knepley     for (i = 0; i < n; ++i) {
105848a46eb9SPierre Jolivet       if (!mglevels[i]->coarseSpace) PetscCall(PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i - 1]->coarseSpace, &mglevels[i]->coarseSpace));
10592b3cbbdaSStefano Zampini       if (i) PetscCall(PCMGAdaptInterpolator_Internal(pc, i, mglevels[i - 1]->smoothu, mglevels[i]->smoothu, mglevels[i - 1]->coarseSpace, mglevels[i]->coarseSpace));
1060f3b08a26SMatthew G. Knepley     }
106148a46eb9SPierre Jolivet     for (i = n - 2; i > -1; --i) PetscCall(PCMGRecomputeLevelOperators_Internal(pc, i));
1062f3b08a26SMatthew G. Knepley   }
1063f3b08a26SMatthew G. Knepley 
10642134b1e4SBarry Smith   if (needRestricts && pc->dm) {
1065caa4e7f2SJed Brown     for (i = n - 2; i >= 0; i--) {
1066ccceb50cSJed Brown       DM  dmfine, dmcoarse;
1067caa4e7f2SJed Brown       Mat Restrict, Inject;
1068caa4e7f2SJed Brown       Vec rscale;
10699566063dSJacob Faibussowitsch       PetscCall(KSPGetDM(mglevels[i + 1]->smoothd, &dmfine));
10709566063dSJacob Faibussowitsch       PetscCall(KSPGetDM(mglevels[i]->smoothd, &dmcoarse));
10719566063dSJacob Faibussowitsch       PetscCall(PCMGGetRestriction(pc, i + 1, &Restrict));
10729566063dSJacob Faibussowitsch       PetscCall(PCMGGetRScale(pc, i + 1, &rscale));
10739566063dSJacob Faibussowitsch       PetscCall(PCMGGetInjection(pc, i + 1, &Inject));
10749566063dSJacob Faibussowitsch       PetscCall(DMRestrict(dmfine, Restrict, rscale, Inject, dmcoarse));
1075caa4e7f2SJed Brown     }
1076caa4e7f2SJed Brown   }
1077f3fbd535SBarry Smith 
1078f3fbd535SBarry Smith   if (!pc->setupcalled) {
107948a46eb9SPierre Jolivet     for (i = 0; i < n; i++) PetscCall(KSPSetFromOptions(mglevels[i]->smoothd));
1080f3fbd535SBarry Smith     for (i = 1; i < n; i++) {
108148a46eb9SPierre Jolivet       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) PetscCall(KSPSetFromOptions(mglevels[i]->smoothu));
108248a46eb9SPierre Jolivet       if (mglevels[i]->cr) PetscCall(KSPSetFromOptions(mglevels[i]->cr));
1083f3fbd535SBarry Smith     }
10843ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
1085f3fbd535SBarry Smith     for (i = 1; i < n; i++) {
10869566063dSJacob Faibussowitsch       PetscCall(PCMGGetInterpolation(pc, i, NULL));
10879566063dSJacob Faibussowitsch       PetscCall(PCMGGetRestriction(pc, i, NULL));
1088f3fbd535SBarry Smith     }
1089f3fbd535SBarry Smith     for (i = 0; i < n - 1; i++) {
1090f3fbd535SBarry Smith       if (!mglevels[i]->b) {
1091f3fbd535SBarry Smith         Vec *vec;
10929566063dSJacob Faibussowitsch         PetscCall(KSPCreateVecs(mglevels[i]->smoothd, 1, &vec, 0, NULL));
10939566063dSJacob Faibussowitsch         PetscCall(PCMGSetRhs(pc, i, *vec));
10949566063dSJacob Faibussowitsch         PetscCall(VecDestroy(vec));
10959566063dSJacob Faibussowitsch         PetscCall(PetscFree(vec));
1096f3fbd535SBarry Smith       }
1097f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
10989566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b, &tvec));
10999566063dSJacob Faibussowitsch         PetscCall(PCMGSetR(pc, i, tvec));
11009566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&tvec));
1101f3fbd535SBarry Smith       }
1102f3fbd535SBarry Smith       if (!mglevels[i]->x) {
11039566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b, &tvec));
11049566063dSJacob Faibussowitsch         PetscCall(PCMGSetX(pc, i, tvec));
11059566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&tvec));
1106f3fbd535SBarry Smith       }
110741b6fd38SMatthew G. Knepley       if (doCR) {
11089566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crx));
11099566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crb));
11109566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(mglevels[i]->crb));
111141b6fd38SMatthew G. Knepley       }
1112f3fbd535SBarry Smith     }
1113f3fbd535SBarry Smith     if (n != 1 && !mglevels[n - 1]->r) {
1114f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
1115f3fbd535SBarry Smith       Vec *vec;
11169566063dSJacob Faibussowitsch       PetscCall(KSPCreateVecs(mglevels[n - 1]->smoothd, 1, &vec, 0, NULL));
11179566063dSJacob Faibussowitsch       PetscCall(PCMGSetR(pc, n - 1, *vec));
11189566063dSJacob Faibussowitsch       PetscCall(VecDestroy(vec));
11199566063dSJacob Faibussowitsch       PetscCall(PetscFree(vec));
1120f3fbd535SBarry Smith     }
112141b6fd38SMatthew G. Knepley     if (doCR) {
11229566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crx));
11239566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crb));
11249566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(mglevels[n - 1]->crb));
112541b6fd38SMatthew G. Knepley     }
1126f3fbd535SBarry Smith   }
1127f3fbd535SBarry Smith 
112898e04984SBarry Smith   if (pc->dm) {
112998e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
113098e04984SBarry Smith     for (i = 0; i < n - 1; i++) {
113198e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
113298e04984SBarry Smith     }
113398e04984SBarry Smith   }
113408549ed4SJed Brown   // We got here (PCSetUp_MG) because the matrix has changed, which means the smoother needs to be set up again (e.g.,
113508549ed4SJed Brown   // new diagonal for Jacobi). Setting it here allows it to be logged under PCSetUp rather than deep inside a PCApply.
113608549ed4SJed Brown   if (mglevels[n - 1]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[n - 1]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
1137f3fbd535SBarry Smith 
1138f3fbd535SBarry Smith   for (i = 1; i < n; i++) {
11392a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1) {
1140f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
11419566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE));
1142f3fbd535SBarry Smith     }
11439566063dSJacob Faibussowitsch     if (mglevels[i]->cr) PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE));
11449566063dSJacob Faibussowitsch     if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
11459566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(mglevels[i]->smoothd));
1146ad540459SPierre Jolivet     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11479566063dSJacob Faibussowitsch     if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
1148d42688cbSBarry Smith     if (!mglevels[i]->residual) {
1149d42688cbSBarry Smith       Mat mat;
11509566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL));
11519566063dSJacob Faibussowitsch       PetscCall(PCMGSetResidual(pc, i, PCMGResidualDefault, mat));
1152d42688cbSBarry Smith     }
1153fcb023d4SJed Brown     if (!mglevels[i]->residualtranspose) {
1154fcb023d4SJed Brown       Mat mat;
11559566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL));
11569566063dSJacob Faibussowitsch       PetscCall(PCMGSetResidualTranspose(pc, i, PCMGResidualTransposeDefault, mat));
1157fcb023d4SJed Brown     }
1158f3fbd535SBarry Smith   }
1159f3fbd535SBarry Smith   for (i = 1; i < n; i++) {
1160f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
1161f3fbd535SBarry Smith       Mat downmat, downpmat;
1162f3fbd535SBarry Smith 
1163f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
11649566063dSJacob Faibussowitsch       PetscCall(KSPGetOperatorsSet(mglevels[i]->smoothu, &opsset, NULL));
1165f3fbd535SBarry Smith       if (!opsset) {
11669566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat));
11679566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[i]->smoothu, downmat, downpmat));
1168f3fbd535SBarry Smith       }
1169f3fbd535SBarry Smith 
11709566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothu, PETSC_TRUE));
11719566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
11729566063dSJacob Faibussowitsch       PetscCall(KSPSetUp(mglevels[i]->smoothu));
1173ad540459SPierre Jolivet       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11749566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
1175f3fbd535SBarry Smith     }
117641b6fd38SMatthew G. Knepley     if (mglevels[i]->cr) {
117741b6fd38SMatthew G. Knepley       Mat downmat, downpmat;
117841b6fd38SMatthew G. Knepley 
117941b6fd38SMatthew G. Knepley       /* check if operators have been set for up, if not use down operators to set them */
11809566063dSJacob Faibussowitsch       PetscCall(KSPGetOperatorsSet(mglevels[i]->cr, &opsset, NULL));
118141b6fd38SMatthew G. Knepley       if (!opsset) {
11829566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat));
11839566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[i]->cr, downmat, downpmat));
118441b6fd38SMatthew G. Knepley       }
118541b6fd38SMatthew G. Knepley 
11869566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE));
11879566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
11889566063dSJacob Faibussowitsch       PetscCall(KSPSetUp(mglevels[i]->cr));
1189ad540459SPierre Jolivet       if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11909566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
119141b6fd38SMatthew G. Knepley     }
1192f3fbd535SBarry Smith   }
1193f3fbd535SBarry Smith 
11949566063dSJacob Faibussowitsch   if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0));
11959566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(mglevels[0]->smoothd));
1196ad540459SPierre Jolivet   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
11979566063dSJacob Faibussowitsch   if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0));
1198f3fbd535SBarry Smith 
1199f3fbd535SBarry Smith     /*
1200f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
1201e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
1202f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
1203f3fbd535SBarry Smith 
1204f3fbd535SBarry Smith    Only support one or the other at the same time.
1205f3fbd535SBarry Smith   */
1206f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
12079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_matlab", &dump, NULL));
1208ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
1209f3fbd535SBarry Smith   dump = PETSC_FALSE;
1210f3fbd535SBarry Smith #endif
12119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_binary", &dump, NULL));
1212ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
1213f3fbd535SBarry Smith 
1214f3fbd535SBarry Smith   if (viewer) {
121548a46eb9SPierre Jolivet     for (i = 1; i < n; i++) PetscCall(MatView(mglevels[i]->restrct, viewer));
1216f3fbd535SBarry Smith     for (i = 0; i < n; i++) {
12179566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(mglevels[i]->smoothd, &pc));
12189566063dSJacob Faibussowitsch       PetscCall(MatView(pc->mat, viewer));
1219f3fbd535SBarry Smith     }
1220f3fbd535SBarry Smith   }
1221f3fbd535SBarry Smith   PetscFunctionReturn(0);
1222f3fbd535SBarry Smith }
1223f3fbd535SBarry Smith 
1224f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
1225f3fbd535SBarry Smith 
1226d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
1227d71ae5a4SJacob Faibussowitsch {
122897d33e41SMatthew G. Knepley   PC_MG *mg = (PC_MG *)pc->data;
122997d33e41SMatthew G. Knepley 
123097d33e41SMatthew G. Knepley   PetscFunctionBegin;
123197d33e41SMatthew G. Knepley   *levels = mg->nlevels;
123297d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
123397d33e41SMatthew G. Knepley }
123497d33e41SMatthew G. Knepley 
12354b9ad928SBarry Smith /*@
1236f1580f4eSBarry Smith    PCMGGetLevels - Gets the number of levels to use with `PCMG`.
12374b9ad928SBarry Smith 
12384b9ad928SBarry Smith    Not Collective
12394b9ad928SBarry Smith 
12404b9ad928SBarry Smith    Input Parameter:
12414b9ad928SBarry Smith .  pc - the preconditioner context
12424b9ad928SBarry Smith 
12434b9ad928SBarry Smith    Output parameter:
12444b9ad928SBarry Smith .  levels - the number of levels
12454b9ad928SBarry Smith 
12464b9ad928SBarry Smith    Level: advanced
12474b9ad928SBarry Smith 
1248f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetLevels()`
12494b9ad928SBarry Smith @*/
1250d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetLevels(PC pc, PetscInt *levels)
1251d71ae5a4SJacob Faibussowitsch {
12524b9ad928SBarry Smith   PetscFunctionBegin;
12530700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12544482741eSBarry Smith   PetscValidIntPointer(levels, 2);
125597d33e41SMatthew G. Knepley   *levels = 0;
1256cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGGetLevels_C", (PC, PetscInt *), (pc, levels));
12574b9ad928SBarry Smith   PetscFunctionReturn(0);
12584b9ad928SBarry Smith }
12594b9ad928SBarry Smith 
12604b9ad928SBarry Smith /*@
1261f1580f4eSBarry Smith    PCMGGetGridComplexity - compute operator and grid complexity of the `PCMG` hierarchy
1262e7d4b4cbSMark Adams 
1263e7d4b4cbSMark Adams    Input Parameter:
1264e7d4b4cbSMark Adams .  pc - the preconditioner context
1265e7d4b4cbSMark Adams 
126690db8557SMark Adams    Output Parameters:
126790db8557SMark Adams +  gc - grid complexity = sum_i(n_i) / n_0
126890db8557SMark Adams -  oc - operator complexity = sum_i(nnz_i) / nnz_0
1269e7d4b4cbSMark Adams 
1270e7d4b4cbSMark Adams    Level: advanced
1271e7d4b4cbSMark Adams 
1272f1580f4eSBarry Smith    Note:
1273f1580f4eSBarry Smith    This is often call the operator complexity in multigrid literature
1274f1580f4eSBarry Smith 
1275f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`
1276e7d4b4cbSMark Adams @*/
1277d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc, PetscReal *oc)
1278d71ae5a4SJacob Faibussowitsch {
1279e7d4b4cbSMark Adams   PC_MG         *mg       = (PC_MG *)pc->data;
1280e7d4b4cbSMark Adams   PC_MG_Levels **mglevels = mg->levels;
1281e7d4b4cbSMark Adams   PetscInt       lev, N;
1282e7d4b4cbSMark Adams   PetscLogDouble nnz0 = 0, sgc = 0, soc = 0, n0 = 0;
1283e7d4b4cbSMark Adams   MatInfo        info;
1284e7d4b4cbSMark Adams 
1285e7d4b4cbSMark Adams   PetscFunctionBegin;
128690db8557SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
128790db8557SMark Adams   if (gc) PetscValidRealPointer(gc, 2);
128890db8557SMark Adams   if (oc) PetscValidRealPointer(oc, 3);
1289e7d4b4cbSMark Adams   if (!pc->setupcalled) {
129090db8557SMark Adams     if (gc) *gc = 0;
129190db8557SMark Adams     if (oc) *oc = 0;
1292e7d4b4cbSMark Adams     PetscFunctionReturn(0);
1293e7d4b4cbSMark Adams   }
1294e7d4b4cbSMark Adams   PetscCheck(mg->nlevels > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MG has no levels");
1295e7d4b4cbSMark Adams   for (lev = 0; lev < mg->nlevels; lev++) {
1296e7d4b4cbSMark Adams     Mat dB;
12979566063dSJacob Faibussowitsch     PetscCall(KSPGetOperators(mglevels[lev]->smoothd, NULL, &dB));
12989566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(dB, MAT_GLOBAL_SUM, &info)); /* global reduction */
12999566063dSJacob Faibussowitsch     PetscCall(MatGetSize(dB, &N, NULL));
1300e7d4b4cbSMark Adams     sgc += N;
1301e7d4b4cbSMark Adams     soc += info.nz_used;
13029371c9d4SSatish Balay     if (lev == mg->nlevels - 1) {
13039371c9d4SSatish Balay       nnz0 = info.nz_used;
13049371c9d4SSatish Balay       n0   = N;
13059371c9d4SSatish Balay     }
1306e7d4b4cbSMark Adams   }
1307*0fdf79fbSJacob Faibussowitsch   PetscCheck(n0 > 0 && gc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number for grid points on finest level is not available");
1308*0fdf79fbSJacob Faibussowitsch   *gc = (PetscReal)(sgc / n0);
130990db8557SMark Adams   if (nnz0 > 0 && oc) *oc = (PetscReal)(soc / nnz0);
1310e7d4b4cbSMark Adams   PetscFunctionReturn(0);
1311e7d4b4cbSMark Adams }
1312e7d4b4cbSMark Adams 
1313e7d4b4cbSMark Adams /*@
131497177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
13154b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
13164b9ad928SBarry Smith 
1317f1580f4eSBarry Smith    Logically Collective on pc
13184b9ad928SBarry Smith 
13194b9ad928SBarry Smith    Input Parameters:
13204b9ad928SBarry Smith +  pc - the preconditioner context
1321f1580f4eSBarry Smith -  form - multigrid form, one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE`
13224b9ad928SBarry Smith 
13234b9ad928SBarry Smith    Options Database Key:
13244b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
13254b9ad928SBarry Smith    additive, full, kaskade
13264b9ad928SBarry Smith 
13274b9ad928SBarry Smith    Level: advanced
13284b9ad928SBarry Smith 
1329f1580f4eSBarry Smith .seealso: `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGGetType()`, `PCMGCycleType`
13304b9ad928SBarry Smith @*/
1331d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetType(PC pc, PCMGType form)
1332d71ae5a4SJacob Faibussowitsch {
1333f3fbd535SBarry Smith   PC_MG *mg = (PC_MG *)pc->data;
13344b9ad928SBarry Smith 
13354b9ad928SBarry Smith   PetscFunctionBegin;
13360700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1337c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc, form, 2);
133831567311SBarry Smith   mg->am = form;
13399dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1340c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
1341c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1342c60c7ad4SBarry Smith }
1343c60c7ad4SBarry Smith 
1344c60c7ad4SBarry Smith /*@
1345f1580f4eSBarry Smith    PCMGGetType - Finds the form of multigrid the `PCMG` is using  multiplicative, additive, full, or the Kaskade algorithm.
1346c60c7ad4SBarry Smith 
1347f1580f4eSBarry Smith    Logically Collective on pc
1348c60c7ad4SBarry Smith 
1349c60c7ad4SBarry Smith    Input Parameter:
1350c60c7ad4SBarry Smith .  pc - the preconditioner context
1351c60c7ad4SBarry Smith 
1352c60c7ad4SBarry Smith    Output Parameter:
1353f1580f4eSBarry Smith .  type - one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE`, `PCMGCycleType`
1354c60c7ad4SBarry Smith 
1355c60c7ad4SBarry Smith    Level: advanced
1356c60c7ad4SBarry Smith 
1357f1580f4eSBarry Smith .seealso: `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGSetType()`
1358c60c7ad4SBarry Smith @*/
1359d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetType(PC pc, PCMGType *type)
1360d71ae5a4SJacob Faibussowitsch {
1361c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG *)pc->data;
1362c60c7ad4SBarry Smith 
1363c60c7ad4SBarry Smith   PetscFunctionBegin;
1364c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1365c60c7ad4SBarry Smith   *type = mg->am;
13664b9ad928SBarry Smith   PetscFunctionReturn(0);
13674b9ad928SBarry Smith }
13684b9ad928SBarry Smith 
13694b9ad928SBarry Smith /*@
1370f1580f4eSBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use `PCMGSetCycleTypeOnLevel()` for more
13714b9ad928SBarry Smith    complicated cycling.
13724b9ad928SBarry Smith 
1373f1580f4eSBarry Smith    Logically Collective on pc
13744b9ad928SBarry Smith 
13754b9ad928SBarry Smith    Input Parameters:
1376c2be2410SBarry Smith +  pc - the multigrid context
1377f1580f4eSBarry Smith -  n - either `PC_MG_CYCLE_V` or `PC_MG_CYCLE_W`
13784b9ad928SBarry Smith 
13794b9ad928SBarry Smith    Options Database Key:
1380c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
13814b9ad928SBarry Smith 
13824b9ad928SBarry Smith    Level: advanced
13834b9ad928SBarry Smith 
1384f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetCycleTypeOnLevel()`, `PCMGType`, `PCMGCycleType`
13854b9ad928SBarry Smith @*/
1386d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetCycleType(PC pc, PCMGCycleType n)
1387d71ae5a4SJacob Faibussowitsch {
1388f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
1389f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
139079416396SBarry Smith   PetscInt       i, levels;
13914b9ad928SBarry Smith 
13924b9ad928SBarry Smith   PetscFunctionBegin;
13930700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
139469fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc, n, 2);
139528b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling");
1396f3fbd535SBarry Smith   levels = mglevels[0]->levels;
13972fa5cd67SKarl Rupp   for (i = 0; i < levels; i++) mglevels[i]->cycles = n;
13984b9ad928SBarry Smith   PetscFunctionReturn(0);
13994b9ad928SBarry Smith }
14004b9ad928SBarry Smith 
14018cc2d5dfSBarry Smith /*@
14028cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
1403f1580f4eSBarry Smith          of multigrid when `PCMGType` is `PC_MG_MULTIPLICATIVE`
14048cc2d5dfSBarry Smith 
1405f1580f4eSBarry Smith    Logically Collective on pc
14068cc2d5dfSBarry Smith 
14078cc2d5dfSBarry Smith    Input Parameters:
14088cc2d5dfSBarry Smith +  pc - the multigrid context
14098cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
14108cc2d5dfSBarry Smith 
14118cc2d5dfSBarry Smith    Options Database Key:
1412147403d9SBarry Smith .  -pc_mg_multiplicative_cycles n - set the number of cycles
14138cc2d5dfSBarry Smith 
14148cc2d5dfSBarry Smith    Level: advanced
14158cc2d5dfSBarry Smith 
1416f1580f4eSBarry Smith    Note:
1417f1580f4eSBarry Smith     This is not associated with setting a v or w cycle, that is set with `PCMGSetCycleType()`
14188cc2d5dfSBarry Smith 
1419f1580f4eSBarry Smith .seealso: `PCMGSetCycleTypeOnLevel()`, `PCMGSetCycleType()`, `PCMGCycleType`, `PCMGType`
14208cc2d5dfSBarry Smith @*/
1421d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGMultiplicativeSetCycles(PC pc, PetscInt n)
1422d71ae5a4SJacob Faibussowitsch {
1423f3fbd535SBarry Smith   PC_MG *mg = (PC_MG *)pc->data;
14248cc2d5dfSBarry Smith 
14258cc2d5dfSBarry Smith   PetscFunctionBegin;
14260700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1427c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
14282a21e185SBarry Smith   mg->cyclesperpcapply = n;
14298cc2d5dfSBarry Smith   PetscFunctionReturn(0);
14308cc2d5dfSBarry Smith }
14318cc2d5dfSBarry Smith 
1432d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetGalerkin_MG(PC pc, PCMGGalerkinType use)
1433d71ae5a4SJacob Faibussowitsch {
1434fb15c04eSBarry Smith   PC_MG *mg = (PC_MG *)pc->data;
1435fb15c04eSBarry Smith 
1436fb15c04eSBarry Smith   PetscFunctionBegin;
14372134b1e4SBarry Smith   mg->galerkin = use;
1438fb15c04eSBarry Smith   PetscFunctionReturn(0);
1439fb15c04eSBarry Smith }
1440fb15c04eSBarry Smith 
1441c2be2410SBarry Smith /*@
144297177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
144382b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1444c2be2410SBarry Smith 
1445f1580f4eSBarry Smith    Logically Collective on pc
1446c2be2410SBarry Smith 
1447c2be2410SBarry Smith    Input Parameters:
1448c91913d3SJed Brown +  pc - the multigrid context
1449f1580f4eSBarry Smith -  use - one of `PC_MG_GALERKIN_BOTH`, `PC_MG_GALERKIN_PMAT`, `PC_MG_GALERKIN_MAT`, or `PC_MG_GALERKIN_NONE`
1450c2be2410SBarry Smith 
1451c2be2410SBarry Smith    Options Database Key:
1452147403d9SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - set the matrices to form via the Galerkin process
1453c2be2410SBarry Smith 
1454c2be2410SBarry Smith    Level: intermediate
1455c2be2410SBarry Smith 
1456f1580f4eSBarry Smith    Note:
1457f1580f4eSBarry Smith    Some codes that use `PCMG` such as `PCGAMG` use Galerkin internally while constructing the hierarchy and thus do not
1458f1580f4eSBarry Smith    use the `PCMG` construction of the coarser grids.
14591c1aac46SBarry Smith 
1460f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetGalerkin()`, `PCMGGalerkinType`
1461c2be2410SBarry Smith @*/
1462d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetGalerkin(PC pc, PCMGGalerkinType use)
1463d71ae5a4SJacob Faibussowitsch {
1464c2be2410SBarry Smith   PetscFunctionBegin;
14650700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1466cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGSetGalerkin_C", (PC, PCMGGalerkinType), (pc, use));
1467c2be2410SBarry Smith   PetscFunctionReturn(0);
1468c2be2410SBarry Smith }
1469c2be2410SBarry Smith 
14703fc8bf9cSBarry Smith /*@
1471f1580f4eSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. A_i-1 = r_i * A_i * p_i
14723fc8bf9cSBarry Smith 
14733fc8bf9cSBarry Smith    Not Collective
14743fc8bf9cSBarry Smith 
14753fc8bf9cSBarry Smith    Input Parameter:
14763fc8bf9cSBarry Smith .  pc - the multigrid context
14773fc8bf9cSBarry Smith 
14783fc8bf9cSBarry Smith    Output Parameter:
1479f1580f4eSBarry 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`
14803fc8bf9cSBarry Smith 
14813fc8bf9cSBarry Smith    Level: intermediate
14823fc8bf9cSBarry Smith 
1483f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetGalerkin()`, `PCMGGalerkinType`
14843fc8bf9cSBarry Smith @*/
1485d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetGalerkin(PC pc, PCMGGalerkinType *galerkin)
1486d71ae5a4SJacob Faibussowitsch {
1487f3fbd535SBarry Smith   PC_MG *mg = (PC_MG *)pc->data;
14883fc8bf9cSBarry Smith 
14893fc8bf9cSBarry Smith   PetscFunctionBegin;
14900700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14912134b1e4SBarry Smith   *galerkin = mg->galerkin;
14923fc8bf9cSBarry Smith   PetscFunctionReturn(0);
14933fc8bf9cSBarry Smith }
14943fc8bf9cSBarry Smith 
1495d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1496d71ae5a4SJacob Faibussowitsch {
1497f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *)pc->data;
1498f3b08a26SMatthew G. Knepley 
1499f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1500f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = adapt;
1501f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1502f3b08a26SMatthew G. Knepley }
1503f3b08a26SMatthew G. Knepley 
1504d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1505d71ae5a4SJacob Faibussowitsch {
1506f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *)pc->data;
1507f3b08a26SMatthew G. Knepley 
1508f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1509f3b08a26SMatthew G. Knepley   *adapt = mg->adaptInterpolation;
1510f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1511f3b08a26SMatthew G. Knepley }
1512f3b08a26SMatthew G. Knepley 
1513d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType ctype)
1514d71ae5a4SJacob Faibussowitsch {
15152b3cbbdaSStefano Zampini   PC_MG *mg = (PC_MG *)pc->data;
15162b3cbbdaSStefano Zampini 
15172b3cbbdaSStefano Zampini   PetscFunctionBegin;
15182b3cbbdaSStefano Zampini   mg->adaptInterpolation = ctype != PCMG_ADAPT_NONE ? PETSC_TRUE : PETSC_FALSE;
15192b3cbbdaSStefano Zampini   mg->coarseSpaceType    = ctype;
15202b3cbbdaSStefano Zampini   PetscFunctionReturn(0);
15212b3cbbdaSStefano Zampini }
15222b3cbbdaSStefano Zampini 
1523d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType *ctype)
1524d71ae5a4SJacob Faibussowitsch {
15252b3cbbdaSStefano Zampini   PC_MG *mg = (PC_MG *)pc->data;
15262b3cbbdaSStefano Zampini 
15272b3cbbdaSStefano Zampini   PetscFunctionBegin;
15282b3cbbdaSStefano Zampini   *ctype = mg->coarseSpaceType;
15292b3cbbdaSStefano Zampini   PetscFunctionReturn(0);
15302b3cbbdaSStefano Zampini }
15312b3cbbdaSStefano Zampini 
1532d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr)
1533d71ae5a4SJacob Faibussowitsch {
153441b6fd38SMatthew G. Knepley   PC_MG *mg = (PC_MG *)pc->data;
153541b6fd38SMatthew G. Knepley 
153641b6fd38SMatthew G. Knepley   PetscFunctionBegin;
153741b6fd38SMatthew G. Knepley   mg->compatibleRelaxation = cr;
153841b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
153941b6fd38SMatthew G. Knepley }
154041b6fd38SMatthew G. Knepley 
1541d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr)
1542d71ae5a4SJacob Faibussowitsch {
154341b6fd38SMatthew G. Knepley   PC_MG *mg = (PC_MG *)pc->data;
154441b6fd38SMatthew G. Knepley 
154541b6fd38SMatthew G. Knepley   PetscFunctionBegin;
154641b6fd38SMatthew G. Knepley   *cr = mg->compatibleRelaxation;
154741b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
154841b6fd38SMatthew G. Knepley }
154941b6fd38SMatthew G. Knepley 
15502b3cbbdaSStefano Zampini /*@C
15512b3cbbdaSStefano Zampini   PCMGSetAdaptCoarseSpaceType - Set the type of adaptive coarse space.
15522b3cbbdaSStefano Zampini 
15532b3cbbdaSStefano 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.
15542b3cbbdaSStefano Zampini 
1555f1580f4eSBarry Smith   Logically Collective on pc
15562b3cbbdaSStefano Zampini 
15572b3cbbdaSStefano Zampini   Input Parameters:
15582b3cbbdaSStefano Zampini + pc    - the multigrid context
15592b3cbbdaSStefano Zampini - ctype - the type of coarse space
15602b3cbbdaSStefano Zampini 
15612b3cbbdaSStefano Zampini   Options Database Keys:
15622b3cbbdaSStefano Zampini + -pc_mg_adapt_interp_n <int>             - The number of modes to use
15632b3cbbdaSStefano Zampini - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: none, polynomial, harmonic, eigenvector, generalized_eigenvector, gdsw
15642b3cbbdaSStefano Zampini 
15652b3cbbdaSStefano Zampini   Level: intermediate
15662b3cbbdaSStefano Zampini 
1567f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGCoarseSpaceType`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()`
15682b3cbbdaSStefano Zampini @*/
1569d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType ctype)
1570d71ae5a4SJacob Faibussowitsch {
15712b3cbbdaSStefano Zampini   PetscFunctionBegin;
15722b3cbbdaSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15732b3cbbdaSStefano Zampini   PetscValidLogicalCollectiveEnum(pc, ctype, 2);
15742b3cbbdaSStefano Zampini   PetscTryMethod(pc, "PCMGSetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType), (pc, ctype));
15752b3cbbdaSStefano Zampini   PetscFunctionReturn(0);
15762b3cbbdaSStefano Zampini }
15772b3cbbdaSStefano Zampini 
15782b3cbbdaSStefano Zampini /*@C
15792b3cbbdaSStefano Zampini    PCMGGetAdaptCoarseSpaceType - Get the type of adaptive coarse space.
15802b3cbbdaSStefano Zampini 
15812b3cbbdaSStefano Zampini    Not Collective
15822b3cbbdaSStefano Zampini 
15832b3cbbdaSStefano Zampini    Input Parameter:
15842b3cbbdaSStefano Zampini .  pc    - the multigrid context
15852b3cbbdaSStefano Zampini 
15862b3cbbdaSStefano Zampini    Output Parameter:
15872b3cbbdaSStefano Zampini .  ctype - the type of coarse space
15882b3cbbdaSStefano Zampini 
15892b3cbbdaSStefano Zampini   Level: intermediate
15902b3cbbdaSStefano Zampini 
1591f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGCoarseSpaceType`, `PCMGSetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()`
15922b3cbbdaSStefano Zampini @*/
1593d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType *ctype)
1594d71ae5a4SJacob Faibussowitsch {
15952b3cbbdaSStefano Zampini   PetscFunctionBegin;
15962b3cbbdaSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15972b3cbbdaSStefano Zampini   PetscValidPointer(ctype, 2);
15982b3cbbdaSStefano Zampini   PetscUseMethod(pc, "PCMGGetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType *), (pc, ctype));
15992b3cbbdaSStefano Zampini   PetscFunctionReturn(0);
16002b3cbbdaSStefano Zampini }
16012b3cbbdaSStefano Zampini 
16022b3cbbdaSStefano Zampini /* MATT: REMOVE? */
1603f3b08a26SMatthew G. Knepley /*@
1604f3b08a26SMatthew 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.
1605f3b08a26SMatthew G. Knepley 
1606f1580f4eSBarry Smith    Logically Collective on pc
1607f3b08a26SMatthew G. Knepley 
1608f3b08a26SMatthew G. Knepley    Input Parameters:
1609f3b08a26SMatthew G. Knepley +  pc    - the multigrid context
1610f3b08a26SMatthew G. Knepley -  adapt - flag for adaptation of the interpolator
1611f3b08a26SMatthew G. Knepley 
1612f3b08a26SMatthew G. Knepley    Options Database Keys:
1613f3b08a26SMatthew G. Knepley +  -pc_mg_adapt_interp                     - Turn on adaptation
1614f3b08a26SMatthew G. Knepley .  -pc_mg_adapt_interp_n <int>             - The number of modes to use, should be divisible by dimension
1615f3b08a26SMatthew G. Knepley -  -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector
1616f3b08a26SMatthew G. Knepley 
1617f3b08a26SMatthew G. Knepley   Level: intermediate
1618f3b08a26SMatthew G. Knepley 
1619f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
1620f3b08a26SMatthew G. Knepley @*/
1621d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1622d71ae5a4SJacob Faibussowitsch {
1623f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1624f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1625cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGSetAdaptInterpolation_C", (PC, PetscBool), (pc, adapt));
1626f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1627f3b08a26SMatthew G. Knepley }
1628f3b08a26SMatthew G. Knepley 
1629f3b08a26SMatthew G. Knepley /*@
1630f1580f4eSBarry Smith   PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh,
1631f1580f4eSBarry Smith   and thus accurately interpolated.
1632f3b08a26SMatthew G. Knepley 
1633f3b08a26SMatthew G. Knepley   Not collective
1634f3b08a26SMatthew G. Knepley 
1635f3b08a26SMatthew G. Knepley   Input Parameter:
1636f3b08a26SMatthew G. Knepley . pc    - the multigrid context
1637f3b08a26SMatthew G. Knepley 
1638f3b08a26SMatthew G. Knepley   Output Parameter:
1639f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator
1640f3b08a26SMatthew G. Knepley 
1641f3b08a26SMatthew G. Knepley   Level: intermediate
1642f3b08a26SMatthew G. Knepley 
1643f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
1644f3b08a26SMatthew G. Knepley @*/
1645d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1646d71ae5a4SJacob Faibussowitsch {
1647f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1648f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1649dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(adapt, 2);
1650cac4c232SBarry Smith   PetscUseMethod(pc, "PCMGGetAdaptInterpolation_C", (PC, PetscBool *), (pc, adapt));
1651f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1652f3b08a26SMatthew G. Knepley }
1653f3b08a26SMatthew G. Knepley 
16544b9ad928SBarry Smith /*@
165541b6fd38SMatthew G. Knepley    PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation.
165641b6fd38SMatthew G. Knepley 
1657f1580f4eSBarry Smith    Logically Collective on pc
165841b6fd38SMatthew G. Knepley 
165941b6fd38SMatthew G. Knepley    Input Parameters:
166041b6fd38SMatthew G. Knepley +  pc - the multigrid context
166141b6fd38SMatthew G. Knepley -  cr - flag for compatible relaxation
166241b6fd38SMatthew G. Knepley 
1663f1580f4eSBarry Smith    Options Database Key:
166441b6fd38SMatthew G. Knepley .  -pc_mg_adapt_cr - Turn on compatible relaxation
166541b6fd38SMatthew G. Knepley 
166641b6fd38SMatthew G. Knepley    Level: intermediate
166741b6fd38SMatthew G. Knepley 
1668f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetAdaptCR()`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
166941b6fd38SMatthew G. Knepley @*/
1670d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr)
1671d71ae5a4SJacob Faibussowitsch {
167241b6fd38SMatthew G. Knepley   PetscFunctionBegin;
167341b6fd38SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1674cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGSetAdaptCR_C", (PC, PetscBool), (pc, cr));
167541b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
167641b6fd38SMatthew G. Knepley }
167741b6fd38SMatthew G. Knepley 
167841b6fd38SMatthew G. Knepley /*@
167941b6fd38SMatthew G. Knepley   PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation.
168041b6fd38SMatthew G. Knepley 
168141b6fd38SMatthew G. Knepley   Not collective
168241b6fd38SMatthew G. Knepley 
168341b6fd38SMatthew G. Knepley   Input Parameter:
168441b6fd38SMatthew G. Knepley . pc    - the multigrid context
168541b6fd38SMatthew G. Knepley 
168641b6fd38SMatthew G. Knepley   Output Parameter:
168741b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion
168841b6fd38SMatthew G. Knepley 
168941b6fd38SMatthew G. Knepley   Level: intermediate
169041b6fd38SMatthew G. Knepley 
1691f1580f4eSBarry Smith .seealso: `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
169241b6fd38SMatthew G. Knepley @*/
1693d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr)
1694d71ae5a4SJacob Faibussowitsch {
169541b6fd38SMatthew G. Knepley   PetscFunctionBegin;
169641b6fd38SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1697dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(cr, 2);
1698cac4c232SBarry Smith   PetscUseMethod(pc, "PCMGGetAdaptCR_C", (PC, PetscBool *), (pc, cr));
169941b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
170041b6fd38SMatthew G. Knepley }
170141b6fd38SMatthew G. Knepley 
170241b6fd38SMatthew G. Knepley /*@
170306792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1704f1580f4eSBarry Smith    on all levels.  Use `PCMGDistinctSmoothUp()` to create separate up and down smoothers if you want different numbers of
1705710315b6SLawrence Mitchell    pre- and post-smoothing steps.
170606792cafSBarry Smith 
1707f1580f4eSBarry Smith    Logically Collective on pc
170806792cafSBarry Smith 
170906792cafSBarry Smith    Input Parameters:
171006792cafSBarry Smith +  mg - the multigrid context
171106792cafSBarry Smith -  n - the number of smoothing steps
171206792cafSBarry Smith 
171306792cafSBarry Smith    Options Database Key:
1714a2b725a8SWilliam Gropp .  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
171506792cafSBarry Smith 
171606792cafSBarry Smith    Level: advanced
171706792cafSBarry Smith 
1718f1580f4eSBarry Smith    Note:
1719f1580f4eSBarry 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.
172006792cafSBarry Smith 
1721f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetDistinctSmoothUp()`
172206792cafSBarry Smith @*/
1723d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetNumberSmooth(PC pc, PetscInt n)
1724d71ae5a4SJacob Faibussowitsch {
172506792cafSBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
172606792cafSBarry Smith   PC_MG_Levels **mglevels = mg->levels;
172706792cafSBarry Smith   PetscInt       i, levels;
172806792cafSBarry Smith 
172906792cafSBarry Smith   PetscFunctionBegin;
173006792cafSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
173106792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
173228b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling");
173306792cafSBarry Smith   levels = mglevels[0]->levels;
173406792cafSBarry Smith 
173506792cafSBarry Smith   for (i = 1; i < levels; i++) {
17369566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[i]->smoothu, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, n));
17379566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, n));
173806792cafSBarry Smith     mg->default_smoothu = n;
173906792cafSBarry Smith     mg->default_smoothd = n;
174006792cafSBarry Smith   }
174106792cafSBarry Smith   PetscFunctionReturn(0);
174206792cafSBarry Smith }
174306792cafSBarry Smith 
1744f442ab6aSBarry Smith /*@
1745f1580f4eSBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate `KSP` from the down (pre) smoother on all levels
1746710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1747f442ab6aSBarry Smith 
1748f1580f4eSBarry Smith    Logically Collective on pc
1749f442ab6aSBarry Smith 
1750f1580f4eSBarry Smith    Input Parameter:
1751f442ab6aSBarry Smith .  pc - the preconditioner context
1752f442ab6aSBarry Smith 
1753f442ab6aSBarry Smith    Options Database Key:
1754147403d9SBarry Smith .  -pc_mg_distinct_smoothup <bool> - use distinct smoothing objects
1755f442ab6aSBarry Smith 
1756f442ab6aSBarry Smith    Level: advanced
1757f442ab6aSBarry Smith 
1758f1580f4eSBarry Smith    Note:
1759f1580f4eSBarry 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.
1760f442ab6aSBarry Smith 
1761f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGSetNumberSmooth()`
1762f442ab6aSBarry Smith @*/
1763d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetDistinctSmoothUp(PC pc)
1764d71ae5a4SJacob Faibussowitsch {
1765f442ab6aSBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
1766f442ab6aSBarry Smith   PC_MG_Levels **mglevels = mg->levels;
1767f442ab6aSBarry Smith   PetscInt       i, levels;
1768f442ab6aSBarry Smith   KSP            subksp;
1769f442ab6aSBarry Smith 
1770f442ab6aSBarry Smith   PetscFunctionBegin;
1771f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
177228b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling");
1773f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1774f442ab6aSBarry Smith 
1775f442ab6aSBarry Smith   for (i = 1; i < levels; i++) {
1776710315b6SLawrence Mitchell     const char *prefix = NULL;
1777f442ab6aSBarry Smith     /* make sure smoother up and down are different */
17789566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmootherUp(pc, i, &subksp));
17799566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(mglevels[i]->smoothd, &prefix));
17809566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(subksp, prefix));
17819566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(subksp, "up_"));
1782f442ab6aSBarry Smith   }
1783f442ab6aSBarry Smith   PetscFunctionReturn(0);
1784f442ab6aSBarry Smith }
1785f442ab6aSBarry Smith 
178607a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1787d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations_MG(PC pc, PetscInt *num_levels, Mat *interpolations[])
1788d71ae5a4SJacob Faibussowitsch {
178907a4832bSFande Kong   PC_MG         *mg       = (PC_MG *)pc->data;
179007a4832bSFande Kong   PC_MG_Levels **mglevels = mg->levels;
179107a4832bSFande Kong   Mat           *mat;
179207a4832bSFande Kong   PetscInt       l;
179307a4832bSFande Kong 
179407a4832bSFande Kong   PetscFunctionBegin;
179528b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
17969566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mg->nlevels, &mat));
179707a4832bSFande Kong   for (l = 1; l < mg->nlevels; l++) {
179807a4832bSFande Kong     mat[l - 1] = mglevels[l]->interpolate;
17999566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat[l - 1]));
180007a4832bSFande Kong   }
180107a4832bSFande Kong   *num_levels     = mg->nlevels;
180207a4832bSFande Kong   *interpolations = mat;
180307a4832bSFande Kong   PetscFunctionReturn(0);
180407a4832bSFande Kong }
180507a4832bSFande Kong 
180607a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1807d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators_MG(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
1808d71ae5a4SJacob Faibussowitsch {
180907a4832bSFande Kong   PC_MG         *mg       = (PC_MG *)pc->data;
181007a4832bSFande Kong   PC_MG_Levels **mglevels = mg->levels;
181107a4832bSFande Kong   PetscInt       l;
181207a4832bSFande Kong   Mat           *mat;
181307a4832bSFande Kong 
181407a4832bSFande Kong   PetscFunctionBegin;
181528b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
18169566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mg->nlevels, &mat));
181707a4832bSFande Kong   for (l = 0; l < mg->nlevels - 1; l++) {
18189566063dSJacob Faibussowitsch     PetscCall(KSPGetOperators(mglevels[l]->smoothd, NULL, &(mat[l])));
18199566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat[l]));
182007a4832bSFande Kong   }
182107a4832bSFande Kong   *num_levels      = mg->nlevels;
182207a4832bSFande Kong   *coarseOperators = mat;
182307a4832bSFande Kong   PetscFunctionReturn(0);
182407a4832bSFande Kong }
182507a4832bSFande Kong 
1826f3b08a26SMatthew G. Knepley /*@C
1827f1580f4eSBarry Smith    PCMGRegisterCoarseSpaceConstructor -  Adds a method to the `PCMG` package for coarse space construction.
1828f3b08a26SMatthew G. Knepley 
1829f3b08a26SMatthew G. Knepley    Not collective
1830f3b08a26SMatthew G. Knepley 
1831f3b08a26SMatthew G. Knepley    Input Parameters:
1832f3b08a26SMatthew G. Knepley +  name     - name of the constructor
1833f3b08a26SMatthew G. Knepley -  function - constructor routine
1834f3b08a26SMatthew G. Knepley 
1835f3b08a26SMatthew G. Knepley    Calling sequence for the routine:
18362b3cbbdaSStefano Zampini $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp)
1837f1580f4eSBarry Smith +  pc        - The PC object
1838f1580f4eSBarry Smith .  l         - The multigrid level, 0 is the coarse level
1839f1580f4eSBarry Smith .  dm        - The DM for this level
1840f1580f4eSBarry Smith .  smooth    - The level smoother
1841f1580f4eSBarry Smith .  Nc        - The size of the coarse space
1842f1580f4eSBarry Smith .  initGuess - Basis for an initial guess for the space
1843f1580f4eSBarry Smith -  coarseSp  - A basis for the computed coarse space
1844f3b08a26SMatthew G. Knepley 
1845f3b08a26SMatthew G. Knepley   Level: advanced
1846f3b08a26SMatthew G. Knepley 
1847f1580f4eSBarry Smith   Developer Note:
1848f1580f4eSBarry Smith   How come this is not used by `PCGAMG`?
1849f1580f4eSBarry Smith 
1850f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGGetCoarseSpaceConstructor()`, `PCRegister()`
1851f3b08a26SMatthew G. Knepley @*/
1852d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *))
1853d71ae5a4SJacob Faibussowitsch {
1854f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
18559566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
18569566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCMGCoarseList, name, function));
1857f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1858f3b08a26SMatthew G. Knepley }
1859f3b08a26SMatthew G. Knepley 
1860f3b08a26SMatthew G. Knepley /*@C
1861f3b08a26SMatthew G. Knepley   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.
1862f3b08a26SMatthew G. Knepley 
1863f3b08a26SMatthew G. Knepley   Not collective
1864f3b08a26SMatthew G. Knepley 
1865f3b08a26SMatthew G. Knepley   Input Parameter:
1866f3b08a26SMatthew G. Knepley . name     - name of the constructor
1867f3b08a26SMatthew G. Knepley 
1868f3b08a26SMatthew G. Knepley   Output Parameter:
1869f3b08a26SMatthew G. Knepley . function - constructor routine
1870f3b08a26SMatthew G. Knepley 
1871f3b08a26SMatthew G. Knepley   Level: advanced
1872f3b08a26SMatthew G. Knepley 
1873f1580f4eSBarry Smith .seealso: `PCMG`, `PCMGRegisterCoarseSpaceConstructor()`, `PCRegister()`
1874f3b08a26SMatthew G. Knepley @*/
1875d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *))
1876d71ae5a4SJacob Faibussowitsch {
1877f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
18789566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PCMGCoarseList, name, function));
1879f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1880f3b08a26SMatthew G. Knepley }
1881f3b08a26SMatthew G. Knepley 
18823b09bd56SBarry Smith /*MC
1883ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
18843b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
18853b09bd56SBarry Smith 
18863b09bd56SBarry Smith    Options Database Keys:
18873b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1888391689abSStefano Zampini .  -pc_mg_cycle_type <v,w> - provide the cycle desired
18898c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
18903b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1891710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
18922134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
18938cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
18948cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1895e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
18968cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
18978cf6037eSBarry Smith                         to the binary output file called binaryoutput
18983b09bd56SBarry Smith 
189995452b02SPatrick Sanan    Notes:
1900f1580f4eSBarry 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
19013b09bd56SBarry Smith 
19028cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
19038cf6037eSBarry Smith 
1904f1580f4eSBarry Smith        When run with `KSPRICHARDSON` the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
190523067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
190623067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
190723067569SBarry Smith        residual is computed at the end of each cycle.
190823067569SBarry Smith 
19093b09bd56SBarry Smith    Level: intermediate
19103b09bd56SBarry Smith 
1911db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCEXOTIC`, `PCGAMG`, `PCML`, `PCHYPRE`
1912db781477SPatrick Sanan           `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `PCMGSetCycleType()`,
1913db781477SPatrick Sanan           `PCMGSetDistinctSmoothUp()`, `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`,
1914db781477SPatrick Sanan           `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`,
1915f1580f4eSBarry Smith           `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`,
1916f1580f4eSBarry Smith           `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
19173b09bd56SBarry Smith M*/
19183b09bd56SBarry Smith 
1919d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1920d71ae5a4SJacob Faibussowitsch {
1921f3fbd535SBarry Smith   PC_MG *mg;
1922f3fbd535SBarry Smith 
19234b9ad928SBarry Smith   PetscFunctionBegin;
19244dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mg));
19253ec1f749SStefano Zampini   pc->data               = mg;
1926f3fbd535SBarry Smith   mg->nlevels            = -1;
192710eca3edSLisandro Dalcin   mg->am                 = PC_MG_MULTIPLICATIVE;
19282134b1e4SBarry Smith   mg->galerkin           = PC_MG_GALERKIN_NONE;
1929f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = PETSC_FALSE;
1930f3b08a26SMatthew G. Knepley   mg->Nc                 = -1;
1931f3b08a26SMatthew G. Knepley   mg->eigenvalue         = -1;
1932f3fbd535SBarry Smith 
193337a44384SMark Adams   pc->useAmat = PETSC_TRUE;
193437a44384SMark Adams 
19354b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
1936fcb023d4SJed Brown   pc->ops->applytranspose = PCApplyTranspose_MG;
193730b0564aSStefano Zampini   pc->ops->matapply       = PCMatApply_MG;
19384b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1939a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
19404b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
19414b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
19424b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1943fb15c04eSBarry Smith 
19449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataRegister(&mg->eigenvalue));
19459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", PCMGSetGalerkin_MG));
19469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", PCMGGetLevels_MG));
19479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", PCMGSetLevels_MG));
19489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", PCGetInterpolations_MG));
19499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", PCGetCoarseOperators_MG));
19509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", PCMGSetAdaptInterpolation_MG));
19519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", PCMGGetAdaptInterpolation_MG));
19529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", PCMGSetAdaptCR_MG));
19539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", PCMGGetAdaptCR_MG));
19542b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", PCMGSetAdaptCoarseSpaceType_MG));
19552b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", PCMGGetAdaptCoarseSpaceType_MG));
19564b9ad928SBarry Smith   PetscFunctionReturn(0);
19574b9ad928SBarry Smith }
1958