xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
34b9ad928SBarry Smith */
4af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/
541b6fd38SMatthew G. Knepley #include <petsc/private/kspimpl.h>
61e25c274SJed Brown #include <petscdm.h>
7391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC, PetscBool *);
84b9ad928SBarry Smith 
9f3b08a26SMatthew G. Knepley /*
10f3b08a26SMatthew G. Knepley    Contains the list of registered coarse space construction routines
11f3b08a26SMatthew G. Knepley */
12f3b08a26SMatthew G. Knepley PetscFunctionList PCMGCoarseList = NULL;
13f3b08a26SMatthew G. Knepley 
14d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGMCycle_Private(PC pc, PC_MG_Levels **mglevelsin, PetscBool transpose, PetscBool matapp, PCRichardsonConvergedReason *reason)
15d71ae5a4SJacob Faibussowitsch {
1631567311SBarry Smith   PC_MG        *mg = (PC_MG *)pc->data;
1731567311SBarry Smith   PC_MG_Levels *mgc, *mglevels = *mglevelsin;
1831567311SBarry Smith   PetscInt      cycles = (mglevels->level == 1) ? 1 : (PetscInt)mglevels->cycles;
194b9ad928SBarry Smith 
204b9ad928SBarry Smith   PetscFunctionBegin;
219566063dSJacob Faibussowitsch   if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve, 0, 0, 0, 0));
22fcb023d4SJed Brown   if (!transpose) {
2330b0564aSStefano Zampini     if (matapp) {
249566063dSJacob Faibussowitsch       PetscCall(KSPMatSolve(mglevels->smoothd, mglevels->B, mglevels->X)); /* pre-smooth */
259566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->smoothd, pc, NULL));
2630b0564aSStefano Zampini     } else {
279566063dSJacob Faibussowitsch       PetscCall(KSPSolve(mglevels->smoothd, mglevels->b, mglevels->x)); /* pre-smooth */
289566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->smoothd, pc, mglevels->x));
2930b0564aSStefano Zampini     }
30fcb023d4SJed Brown   } else {
3128b400f6SJacob Faibussowitsch     PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported");
329566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(mglevels->smoothu, mglevels->b, mglevels->x)); /* transpose of post-smooth */
339566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(mglevels->smoothu, pc, mglevels->x));
34fcb023d4SJed Brown   }
359566063dSJacob Faibussowitsch   if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve, 0, 0, 0, 0));
3631567311SBarry Smith   if (mglevels->level) { /* not the coarsest grid */
379566063dSJacob Faibussowitsch     if (mglevels->eventresidual) PetscCall(PetscLogEventBegin(mglevels->eventresidual, 0, 0, 0, 0));
3848a46eb9SPierre Jolivet     if (matapp && !mglevels->R) PetscCall(MatDuplicate(mglevels->B, MAT_DO_NOT_COPY_VALUES, &mglevels->R));
39fcb023d4SJed Brown     if (!transpose) {
409566063dSJacob Faibussowitsch       if (matapp) PetscCall((*mglevels->matresidual)(mglevels->A, mglevels->B, mglevels->X, mglevels->R));
419566063dSJacob Faibussowitsch       else PetscCall((*mglevels->residual)(mglevels->A, mglevels->b, mglevels->x, mglevels->r));
42fcb023d4SJed Brown     } else {
439566063dSJacob Faibussowitsch       if (matapp) PetscCall((*mglevels->matresidualtranspose)(mglevels->A, mglevels->B, mglevels->X, mglevels->R));
449566063dSJacob Faibussowitsch       else PetscCall((*mglevels->residualtranspose)(mglevels->A, mglevels->b, mglevels->x, mglevels->r));
45fcb023d4SJed Brown     }
469566063dSJacob Faibussowitsch     if (mglevels->eventresidual) PetscCall(PetscLogEventEnd(mglevels->eventresidual, 0, 0, 0, 0));
474b9ad928SBarry Smith 
484b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
4931567311SBarry Smith     if (mglevels->level == mglevels->levels - 1 && mg->ttol && reason) {
504b9ad928SBarry Smith       PetscReal rnorm;
519566063dSJacob Faibussowitsch       PetscCall(VecNorm(mglevels->r, NORM_2, &rnorm));
524b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
5370441072SBarry Smith         if (rnorm < mg->abstol) {
544d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
559566063dSJacob Faibussowitsch           PetscCall(PetscInfo(pc, "Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n", (double)rnorm, (double)mg->abstol));
564b9ad928SBarry Smith         } else {
574d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
589566063dSJacob Faibussowitsch           PetscCall(PetscInfo(pc, "Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n", (double)rnorm, (double)mg->ttol));
594b9ad928SBarry Smith         }
603ba16761SJacob Faibussowitsch         PetscFunctionReturn(PETSC_SUCCESS);
614b9ad928SBarry Smith       }
624b9ad928SBarry Smith     }
634b9ad928SBarry Smith 
6431567311SBarry Smith     mgc = *(mglevelsin - 1);
659566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict, 0, 0, 0, 0));
66fcb023d4SJed Brown     if (!transpose) {
679566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatRestrict(mglevels->restrct, mglevels->R, &mgc->B));
689566063dSJacob Faibussowitsch       else PetscCall(MatRestrict(mglevels->restrct, mglevels->r, mgc->b));
69fcb023d4SJed Brown     } else {
709566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatRestrict(mglevels->interpolate, mglevels->R, &mgc->B));
719566063dSJacob Faibussowitsch       else PetscCall(MatRestrict(mglevels->interpolate, mglevels->r, mgc->b));
72fcb023d4SJed Brown     }
739566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict, 0, 0, 0, 0));
7430b0564aSStefano Zampini     if (matapp) {
7530b0564aSStefano Zampini       if (!mgc->X) {
769566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(mgc->B, MAT_DO_NOT_COPY_VALUES, &mgc->X));
7730b0564aSStefano Zampini       } else {
789566063dSJacob Faibussowitsch         PetscCall(MatZeroEntries(mgc->X));
7930b0564aSStefano Zampini       }
8030b0564aSStefano Zampini     } else {
819566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(mgc->x));
8230b0564aSStefano Zampini     }
8348a46eb9SPierre Jolivet     while (cycles--) PetscCall(PCMGMCycle_Private(pc, mglevelsin - 1, transpose, matapp, reason));
849566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict, 0, 0, 0, 0));
85fcb023d4SJed Brown     if (!transpose) {
869566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatInterpolateAdd(mglevels->interpolate, mgc->X, mglevels->X, &mglevels->X));
879566063dSJacob Faibussowitsch       else PetscCall(MatInterpolateAdd(mglevels->interpolate, mgc->x, mglevels->x, mglevels->x));
88fcb023d4SJed Brown     } else {
899566063dSJacob Faibussowitsch       PetscCall(MatInterpolateAdd(mglevels->restrct, mgc->x, mglevels->x, mglevels->x));
90fcb023d4SJed Brown     }
919566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict, 0, 0, 0, 0));
929566063dSJacob Faibussowitsch     if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve, 0, 0, 0, 0));
93fcb023d4SJed Brown     if (!transpose) {
9430b0564aSStefano Zampini       if (matapp) {
959566063dSJacob Faibussowitsch         PetscCall(KSPMatSolve(mglevels->smoothu, mglevels->B, mglevels->X)); /* post smooth */
969566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(mglevels->smoothu, pc, NULL));
9730b0564aSStefano Zampini       } else {
989566063dSJacob Faibussowitsch         PetscCall(KSPSolve(mglevels->smoothu, mglevels->b, mglevels->x)); /* post smooth */
999566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(mglevels->smoothu, pc, mglevels->x));
10030b0564aSStefano Zampini       }
101fcb023d4SJed Brown     } else {
10228b400f6SJacob Faibussowitsch       PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported");
1039566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(mglevels->smoothd, mglevels->b, mglevels->x)); /* post smooth */
1049566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->smoothd, pc, mglevels->x));
105fcb023d4SJed Brown     }
10641b6fd38SMatthew G. Knepley     if (mglevels->cr) {
10728b400f6SJacob Faibussowitsch       PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported");
10841b6fd38SMatthew G. Knepley       /* TODO Turn on copy and turn off noisy if we have an exact solution
1099566063dSJacob Faibussowitsch       PetscCall(VecCopy(mglevels->x, mglevels->crx));
1109566063dSJacob Faibussowitsch       PetscCall(VecCopy(mglevels->b, mglevels->crb)); */
1119566063dSJacob Faibussowitsch       PetscCall(KSPSetNoisy_Private(mglevels->crx));
1129566063dSJacob Faibussowitsch       PetscCall(KSPSolve(mglevels->cr, mglevels->crb, mglevels->crx)); /* compatible relaxation */
1139566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->cr, pc, mglevels->crx));
11441b6fd38SMatthew G. Knepley     }
1159566063dSJacob Faibussowitsch     if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve, 0, 0, 0, 0));
1164b9ad928SBarry Smith   }
1173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1184b9ad928SBarry Smith }
1194b9ad928SBarry Smith 
120d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyRichardson_MG(PC pc, Vec b, Vec x, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool zeroguess, PetscInt *outits, PCRichardsonConvergedReason *reason)
121d71ae5a4SJacob Faibussowitsch {
122f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
123f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
124391689abSStefano Zampini   PC             tpc;
125391689abSStefano Zampini   PetscBool      changeu, changed;
126f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels, i;
1274b9ad928SBarry Smith 
1284b9ad928SBarry Smith   PetscFunctionBegin;
129a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
130a762d673SBarry Smith   for (i = 0; i < levels; i++) {
131a762d673SBarry Smith     if (!mglevels[i]->A) {
1329566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL));
1339566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A));
134a762d673SBarry Smith     }
135a762d673SBarry Smith   }
136391689abSStefano Zampini 
1379566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc));
1389566063dSJacob Faibussowitsch   PetscCall(PCPreSolveChangeRHS(tpc, &changed));
1399566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc));
1409566063dSJacob Faibussowitsch   PetscCall(PCPreSolveChangeRHS(tpc, &changeu));
141391689abSStefano Zampini   if (!changed && !changeu) {
1429566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[levels - 1]->b));
143f3fbd535SBarry Smith     mglevels[levels - 1]->b = b;
144391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
145391689abSStefano Zampini     if (!mglevels[levels - 1]->b) {
146391689abSStefano Zampini       Vec *vec;
147391689abSStefano Zampini 
1489566063dSJacob Faibussowitsch       PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL));
149391689abSStefano Zampini       mglevels[levels - 1]->b = *vec;
1509566063dSJacob Faibussowitsch       PetscCall(PetscFree(vec));
151391689abSStefano Zampini     }
1529566063dSJacob Faibussowitsch     PetscCall(VecCopy(b, mglevels[levels - 1]->b));
153391689abSStefano Zampini   }
154f3fbd535SBarry Smith   mglevels[levels - 1]->x = x;
1554b9ad928SBarry Smith 
15631567311SBarry Smith   mg->rtol   = rtol;
15731567311SBarry Smith   mg->abstol = abstol;
15831567311SBarry Smith   mg->dtol   = dtol;
1594b9ad928SBarry Smith   if (rtol) {
1604b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
1614b9ad928SBarry Smith     PetscReal rnorm;
1627319c654SBarry Smith     if (zeroguess) {
1639566063dSJacob Faibussowitsch       PetscCall(VecNorm(b, NORM_2, &rnorm));
1647319c654SBarry Smith     } else {
1659566063dSJacob Faibussowitsch       PetscCall((*mglevels[levels - 1]->residual)(mglevels[levels - 1]->A, b, x, w));
1669566063dSJacob Faibussowitsch       PetscCall(VecNorm(w, NORM_2, &rnorm));
1677319c654SBarry Smith     }
16831567311SBarry Smith     mg->ttol = PetscMax(rtol * rnorm, abstol);
1692fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
1702fa5cd67SKarl Rupp   else mg->ttol = 0.0;
1714b9ad928SBarry Smith 
1724d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
173336babb1SJed Brown      stop prematurely due to small residual */
1744d0a8057SBarry Smith   for (i = 1; i < levels; i++) {
1759566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[i]->smoothu, 0, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
176f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
17723067569SBarry Smith       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
1789566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE));
1799566063dSJacob Faibussowitsch       PetscCall(KSPSetTolerances(mglevels[i]->smoothd, 0, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
1804b9ad928SBarry Smith     }
1814d0a8057SBarry Smith   }
1824d0a8057SBarry Smith 
1834d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1844d0a8057SBarry Smith   for (i = 0; i < its; i++) {
1859566063dSJacob Faibussowitsch     PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, PETSC_FALSE, PETSC_FALSE, reason));
1864d0a8057SBarry Smith     if (*reason) break;
1874d0a8057SBarry Smith   }
1884d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1894d0a8057SBarry Smith   *outits = i;
190391689abSStefano Zampini   if (!changed && !changeu) mglevels[levels - 1]->b = NULL;
1913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1924b9ad928SBarry Smith }
1934b9ad928SBarry Smith 
194d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset_MG(PC pc)
195d71ae5a4SJacob Faibussowitsch {
1963aeaf226SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
1973aeaf226SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
1982b3cbbdaSStefano Zampini   PetscInt       i, n;
1993aeaf226SBarry Smith 
2003aeaf226SBarry Smith   PetscFunctionBegin;
2013aeaf226SBarry Smith   if (mglevels) {
2023aeaf226SBarry Smith     n = mglevels[0]->levels;
2033aeaf226SBarry Smith     for (i = 0; i < n - 1; i++) {
2049566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i + 1]->r));
2059566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->b));
2069566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->x));
2079566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i + 1]->R));
2089566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i]->B));
2099566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i]->X));
2109566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->crx));
2119566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->crb));
2129566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i + 1]->restrct));
2139566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i + 1]->interpolate));
2149566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i + 1]->inject));
2159566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i + 1]->rscale));
2163aeaf226SBarry Smith     }
2179566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[n - 1]->crx));
2189566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[n - 1]->crb));
219391689abSStefano Zampini     /* this is not null only if the smoother on the finest level
220391689abSStefano Zampini        changes the rhs during PreSolve */
2219566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[n - 1]->b));
2229566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&mglevels[n - 1]->B));
2233aeaf226SBarry Smith 
2243aeaf226SBarry Smith     for (i = 0; i < n; i++) {
2252b3cbbdaSStefano Zampini       PetscCall(MatDestroy(&mglevels[i]->coarseSpace));
2269566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i]->A));
22748a46eb9SPierre Jolivet       if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPReset(mglevels[i]->smoothd));
2289566063dSJacob Faibussowitsch       PetscCall(KSPReset(mglevels[i]->smoothu));
2299566063dSJacob Faibussowitsch       if (mglevels[i]->cr) PetscCall(KSPReset(mglevels[i]->cr));
2303aeaf226SBarry Smith     }
231f3b08a26SMatthew G. Knepley     mg->Nc = 0;
2323aeaf226SBarry Smith   }
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2343aeaf226SBarry Smith }
2353aeaf226SBarry Smith 
23641b6fd38SMatthew G. Knepley /* Implementing CR
23741b6fd38SMatthew G. Knepley 
23841b6fd38SMatthew G. Knepley We only want to make corrections that ``do not change'' the coarse solution. What we mean by not changing is that if I prolong my coarse solution to the fine grid and then inject that fine solution back to the coarse grid, I get the same answer. Injection is what Brannick calls R. We want the complementary projector to Inj, which we will call S, after Brannick, so that Inj S = 0. Now the orthogonal projector onto the range of Inj^T is
23941b6fd38SMatthew G. Knepley 
24041b6fd38SMatthew G. Knepley   Inj^T (Inj Inj^T)^{-1} Inj
24141b6fd38SMatthew G. Knepley 
24241b6fd38SMatthew G. Knepley and if Inj is a VecScatter, as it is now in PETSc, we have
24341b6fd38SMatthew G. Knepley 
24441b6fd38SMatthew G. Knepley   Inj^T Inj
24541b6fd38SMatthew G. Knepley 
24641b6fd38SMatthew G. Knepley and
24741b6fd38SMatthew G. Knepley 
24841b6fd38SMatthew G. Knepley   S = I - Inj^T Inj
24941b6fd38SMatthew G. Knepley 
25041b6fd38SMatthew G. Knepley since
25141b6fd38SMatthew G. Knepley 
25241b6fd38SMatthew G. Knepley   Inj S = Inj - (Inj Inj^T) Inj = 0.
25341b6fd38SMatthew G. Knepley 
25441b6fd38SMatthew G. Knepley Brannick suggests
25541b6fd38SMatthew G. Knepley 
25641b6fd38SMatthew G. Knepley   A \to S^T A S  \qquad\mathrm{and}\qquad M \to S^T M S
25741b6fd38SMatthew G. Knepley 
25841b6fd38SMatthew G. Knepley but I do not think his :math:`S^T S = I` is correct. Our S is an orthogonal projector, so :math:`S^T S = S^2 = S`. We will use
25941b6fd38SMatthew G. Knepley 
26041b6fd38SMatthew G. Knepley   M^{-1} A \to S M^{-1} A S
26141b6fd38SMatthew G. Knepley 
26241b6fd38SMatthew G. Knepley In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left.
26341b6fd38SMatthew G. Knepley 
26441b6fd38SMatthew G. Knepley   Check: || Inj P - I ||_F < tol
26541b6fd38SMatthew G. Knepley   Check: In general, Inj Inj^T = I
26641b6fd38SMatthew G. Knepley */
26741b6fd38SMatthew G. Knepley 
26841b6fd38SMatthew G. Knepley typedef struct {
26941b6fd38SMatthew G. Knepley   PC       mg;  /* The PCMG object */
27041b6fd38SMatthew G. Knepley   PetscInt l;   /* The multigrid level for this solver */
27141b6fd38SMatthew G. Knepley   Mat      Inj; /* The injection matrix */
27241b6fd38SMatthew G. Knepley   Mat      S;   /* I - Inj^T Inj */
27341b6fd38SMatthew G. Knepley } CRContext;
27441b6fd38SMatthew G. Knepley 
275d71ae5a4SJacob Faibussowitsch static PetscErrorCode CRSetup_Private(PC pc)
276d71ae5a4SJacob Faibussowitsch {
27741b6fd38SMatthew G. Knepley   CRContext *ctx;
27841b6fd38SMatthew G. Knepley   Mat        It;
27941b6fd38SMatthew G. Knepley 
28041b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
2819566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &ctx));
2829566063dSJacob Faibussowitsch   PetscCall(PCMGGetInjection(ctx->mg, ctx->l, &It));
28328b400f6SJacob Faibussowitsch   PetscCheck(It, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG");
2849566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(It, &ctx->Inj));
2859566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(ctx->Inj, &ctx->S));
2869566063dSJacob Faibussowitsch   PetscCall(MatScale(ctx->S, -1.0));
2879566063dSJacob Faibussowitsch   PetscCall(MatShift(ctx->S, 1.0));
2883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28941b6fd38SMatthew G. Knepley }
29041b6fd38SMatthew G. Knepley 
291d71ae5a4SJacob Faibussowitsch static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y)
292d71ae5a4SJacob Faibussowitsch {
29341b6fd38SMatthew G. Knepley   CRContext *ctx;
29441b6fd38SMatthew G. Knepley 
29541b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
2969566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &ctx));
2979566063dSJacob Faibussowitsch   PetscCall(MatMult(ctx->S, x, y));
2983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
29941b6fd38SMatthew G. Knepley }
30041b6fd38SMatthew G. Knepley 
301d71ae5a4SJacob Faibussowitsch static PetscErrorCode CRDestroy_Private(PC pc)
302d71ae5a4SJacob Faibussowitsch {
30341b6fd38SMatthew G. Knepley   CRContext *ctx;
30441b6fd38SMatthew G. Knepley 
30541b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
3069566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &ctx));
3079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->Inj));
3089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->S));
3099566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
3109566063dSJacob Faibussowitsch   PetscCall(PCShellSetContext(pc, NULL));
3113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
31241b6fd38SMatthew G. Knepley }
31341b6fd38SMatthew G. Knepley 
314d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr)
315d71ae5a4SJacob Faibussowitsch {
31641b6fd38SMatthew G. Knepley   CRContext *ctx;
31741b6fd38SMatthew G. Knepley 
31841b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
3199566063dSJacob Faibussowitsch   PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), cr));
3209566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*cr, "S (complementary projector to injection)"));
3219566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(1, &ctx));
32241b6fd38SMatthew G. Knepley   ctx->mg = pc;
32341b6fd38SMatthew G. Knepley   ctx->l  = l;
3249566063dSJacob Faibussowitsch   PetscCall(PCSetType(*cr, PCSHELL));
3259566063dSJacob Faibussowitsch   PetscCall(PCShellSetContext(*cr, ctx));
3269566063dSJacob Faibussowitsch   PetscCall(PCShellSetApply(*cr, CRApply_Private));
3279566063dSJacob Faibussowitsch   PetscCall(PCShellSetSetUp(*cr, CRSetup_Private));
3289566063dSJacob Faibussowitsch   PetscCall(PCShellSetDestroy(*cr, CRDestroy_Private));
3293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33041b6fd38SMatthew G. Knepley }
33141b6fd38SMatthew G. Knepley 
3328f4fb22eSMark Adams PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions, const char[], const char[], const char *[], const char *[], PetscBool *);
3338f4fb22eSMark Adams 
334d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetLevels_MG(PC pc, PetscInt levels, MPI_Comm *comms)
335d71ae5a4SJacob Faibussowitsch {
336f3fbd535SBarry Smith   PC_MG         *mg = (PC_MG *)pc->data;
337ce94432eSBarry Smith   MPI_Comm       comm;
3383aeaf226SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
33910eca3edSLisandro Dalcin   PCMGType       mgtype   = mg->am;
34010167fecSBarry Smith   PetscInt       mgctype  = (PetscInt)PC_MG_CYCLE_V;
341f3fbd535SBarry Smith   PetscInt       i;
342f3fbd535SBarry Smith   PetscMPIInt    size;
343f3fbd535SBarry Smith   const char    *prefix;
344f3fbd535SBarry Smith   PC             ipc;
3453aeaf226SBarry Smith   PetscInt       n;
3464b9ad928SBarry Smith 
3474b9ad928SBarry Smith   PetscFunctionBegin;
3480700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
349548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc, levels, 2);
3503ba16761SJacob Faibussowitsch   if (mg->nlevels == levels) PetscFunctionReturn(PETSC_SUCCESS);
3519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
3523aeaf226SBarry Smith   if (mglevels) {
35310eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
3543aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
3559566063dSJacob Faibussowitsch     PetscCall(PCReset_MG(pc));
3563aeaf226SBarry Smith     n = mglevels[0]->levels;
3573aeaf226SBarry Smith     for (i = 0; i < n; i++) {
35848a46eb9SPierre Jolivet       if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPDestroy(&mglevels[i]->smoothd));
3599566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&mglevels[i]->smoothu));
3609566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&mglevels[i]->cr));
3619566063dSJacob Faibussowitsch       PetscCall(PetscFree(mglevels[i]));
3623aeaf226SBarry Smith     }
3639566063dSJacob Faibussowitsch     PetscCall(PetscFree(mg->levels));
3643aeaf226SBarry Smith   }
365f3fbd535SBarry Smith 
366f3fbd535SBarry Smith   mg->nlevels = levels;
367f3fbd535SBarry Smith 
3689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(levels, &mglevels));
369f3fbd535SBarry Smith 
3709566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &prefix));
371f3fbd535SBarry Smith 
372a9db87a2SMatthew G Knepley   mg->stageApply = 0;
373f3fbd535SBarry Smith   for (i = 0; i < levels; i++) {
3744dfa11a4SJacob Faibussowitsch     PetscCall(PetscNew(&mglevels[i]));
3752fa5cd67SKarl Rupp 
376f3fbd535SBarry Smith     mglevels[i]->level               = i;
377f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
37810eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
379336babb1SJed Brown     mg->default_smoothu              = 2;
380336babb1SJed Brown     mg->default_smoothd              = 2;
38163e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
38263e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
38363e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
38463e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
385f3fbd535SBarry Smith 
386f3fbd535SBarry Smith     if (comms) comm = comms[i];
387d5a8f7e6SBarry Smith     if (comm != MPI_COMM_NULL) {
3889566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm, &mglevels[i]->smoothd));
3893821be0aSBarry Smith       PetscCall(KSPSetNestLevel(mglevels[i]->smoothd, pc->kspnestlevel));
3909566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->smoothd, pc->erroriffailure));
3919566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd, (PetscObject)pc, levels - i));
3929566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd, prefix));
3939566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level));
3948f4fb22eSMark Adams       if (i == 0 && levels > 1) { // coarse grid
3959566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(mglevels[0]->smoothd, "mg_coarse_"));
396f3fbd535SBarry Smith 
3979dbfc187SHong Zhang         /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
3989566063dSJacob Faibussowitsch         PetscCall(KSPSetType(mglevels[0]->smoothd, KSPPREONLY));
3999566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(mglevels[0]->smoothd, &ipc));
4009566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(comm, &size));
401f3fbd535SBarry Smith         if (size > 1) {
4029566063dSJacob Faibussowitsch           PetscCall(PCSetType(ipc, PCREDUNDANT));
403f3fbd535SBarry Smith         } else {
4049566063dSJacob Faibussowitsch           PetscCall(PCSetType(ipc, PCLU));
405f3fbd535SBarry Smith         }
4069566063dSJacob Faibussowitsch         PetscCall(PCFactorSetShiftType(ipc, MAT_SHIFT_INBLOCKS));
4078f4fb22eSMark Adams       } else {
4088f4fb22eSMark Adams         char tprefix[128];
4098f4fb22eSMark Adams 
4108f4fb22eSMark Adams         PetscCall(KSPSetType(mglevels[i]->smoothd, KSPCHEBYSHEV));
4118f4fb22eSMark Adams         PetscCall(KSPSetConvergenceTest(mglevels[i]->smoothd, KSPConvergedSkip, NULL, NULL));
4128f4fb22eSMark Adams         PetscCall(KSPSetNormType(mglevels[i]->smoothd, KSP_NORM_NONE));
4138f4fb22eSMark Adams         PetscCall(KSPGetPC(mglevels[i]->smoothd, &ipc));
4148f4fb22eSMark Adams         PetscCall(PCSetType(ipc, PCSOR));
4158f4fb22eSMark Adams         PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd));
4168f4fb22eSMark Adams 
4178f4fb22eSMark Adams         if (i == levels - 1 && levels > 1) { // replace 'mg_finegrid_' with 'mg_levels_X_'
4188f4fb22eSMark Adams           PetscBool set;
4198f4fb22eSMark Adams           PetscCall(PetscOptionsFindPairPrefix_Private(((PetscObject)mglevels[i]->smoothd)->options, ((PetscObject)mglevels[i]->smoothd)->prefix, "-mg_fine_", NULL, NULL, &set));
4208f4fb22eSMark Adams           if (set) {
4218f4fb22eSMark Adams             if (prefix) PetscCall(PetscSNPrintf(tprefix, 128, "%smg_fine_", prefix));
4228f4fb22eSMark Adams             else PetscCall(PetscSNPrintf(tprefix, 128, "mg_fine_"));
4238f4fb22eSMark Adams             PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd, tprefix));
4248f4fb22eSMark Adams           } else {
4258f4fb22eSMark Adams             PetscCall(PetscSNPrintf(tprefix, 128, "mg_levels_%d_", (int)i));
4268f4fb22eSMark Adams             PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd, tprefix));
4278f4fb22eSMark Adams           }
4288f4fb22eSMark Adams         } else {
4298f4fb22eSMark Adams           PetscCall(PetscSNPrintf(tprefix, 128, "mg_levels_%d_", (int)i));
4308f4fb22eSMark Adams           PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd, tprefix));
4318f4fb22eSMark Adams         }
432f3fbd535SBarry Smith       }
433d5a8f7e6SBarry Smith     }
434f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
43531567311SBarry Smith     mg->rtol             = 0.0;
43631567311SBarry Smith     mg->abstol           = 0.0;
43731567311SBarry Smith     mg->dtol             = 0.0;
43831567311SBarry Smith     mg->ttol             = 0.0;
43931567311SBarry Smith     mg->cyclesperpcapply = 1;
440f3fbd535SBarry Smith   }
441f3fbd535SBarry Smith   mg->levels = mglevels;
4429566063dSJacob Faibussowitsch   PetscCall(PCMGSetType(pc, mgtype));
4433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4444b9ad928SBarry Smith }
4454b9ad928SBarry Smith 
44697d33e41SMatthew G. Knepley /*@C
447f1580f4eSBarry Smith   PCMGSetLevels - Sets the number of levels to use with `PCMG`.
448f1580f4eSBarry Smith   Must be called before any other `PCMG` routine.
44997d33e41SMatthew G. Knepley 
450c3339decSBarry Smith   Logically Collective
45197d33e41SMatthew G. Knepley 
45297d33e41SMatthew G. Knepley   Input Parameters:
45397d33e41SMatthew G. Knepley + pc     - the preconditioner context
45497d33e41SMatthew G. Knepley . levels - the number of levels
45597d33e41SMatthew G. Knepley - comms  - optional communicators for each level; this is to allow solving the coarser problems
456d5a8f7e6SBarry Smith            on smaller sets of processes. For processes that are not included in the computation
45720f4b53cSBarry Smith            you must pass `MPI_COMM_NULL`. Use comms = `NULL` to specify that all processes
45805552d4bSJunchao Zhang            should participate in each level of problem.
45997d33e41SMatthew G. Knepley 
46097d33e41SMatthew G. Knepley   Level: intermediate
46197d33e41SMatthew G. Knepley 
46297d33e41SMatthew G. Knepley   Notes:
46320f4b53cSBarry Smith   If the number of levels is one then the multigrid uses the `-mg_levels` prefix
4648f4fb22eSMark Adams   for setting the level options rather than the `-mg_coarse` or `-mg_fine` prefix.
46597d33e41SMatthew G. Knepley 
466d5a8f7e6SBarry Smith   You can free the information in comms after this routine is called.
467d5a8f7e6SBarry Smith 
468f1580f4eSBarry Smith   The array of MPI communicators must contain `MPI_COMM_NULL` for those ranks that at each level
469d5a8f7e6SBarry Smith   are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on
470d5a8f7e6SBarry Smith   the two levels, rank 0 in the original communicator will pass in an array of 2 communicators
471d5a8f7e6SBarry Smith   of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators
472f1580f4eSBarry Smith   the first of size 2 and the second of value `MPI_COMM_NULL` since the rank 1 does not participate
473d5a8f7e6SBarry Smith   in the coarse grid solve.
474d5a8f7e6SBarry Smith 
475f1580f4eSBarry Smith   Since each coarser level may have a new `MPI_Comm` with fewer ranks than the previous, one
476d5a8f7e6SBarry Smith   must take special care in providing the restriction and interpolation operation. We recommend
477d5a8f7e6SBarry Smith   providing these as two step operations; first perform a standard restriction or interpolation on
478d5a8f7e6SBarry Smith   the full number of ranks for that level and then use an MPI call to copy the resulting vector
47905552d4bSJunchao Zhang   array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, note in both
480d5a8f7e6SBarry Smith   cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and
48120f4b53cSBarry Smith   receives or `MPI_AlltoAllv()` could be used to do the reshuffling of the vector entries.
482d5a8f7e6SBarry Smith 
483feefa0e1SJacob Faibussowitsch   Fortran Notes:
48420f4b53cSBarry Smith   Use comms = `PETSC_NULL_MPI_COMM` as the equivalent of `NULL` in the C interface. Note `PETSC_NULL_MPI_COMM`
485f1580f4eSBarry Smith   is not `MPI_COMM_NULL`. It is more like `PETSC_NULL_INTEGER`, `PETSC_NULL_REAL` etc.
486d5a8f7e6SBarry Smith 
487562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGSetType()`, `PCMGGetLevels()`
48897d33e41SMatthew G. Knepley @*/
489d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetLevels(PC pc, PetscInt levels, MPI_Comm *comms)
490d71ae5a4SJacob Faibussowitsch {
49197d33e41SMatthew G. Knepley   PetscFunctionBegin;
49297d33e41SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
4934f572ea9SToby Isaac   if (comms) PetscAssertPointer(comms, 3);
494cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGSetLevels_C", (PC, PetscInt, MPI_Comm *), (pc, levels, comms));
4953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
49697d33e41SMatthew G. Knepley }
49797d33e41SMatthew G. Knepley 
498d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy_MG(PC pc)
499d71ae5a4SJacob Faibussowitsch {
500a06653b4SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
501a06653b4SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
502a06653b4SBarry Smith   PetscInt       i, n;
503c07bf074SBarry Smith 
504c07bf074SBarry Smith   PetscFunctionBegin;
5059566063dSJacob Faibussowitsch   PetscCall(PCReset_MG(pc));
506a06653b4SBarry Smith   if (mglevels) {
507a06653b4SBarry Smith     n = mglevels[0]->levels;
508a06653b4SBarry Smith     for (i = 0; i < n; i++) {
50948a46eb9SPierre Jolivet       if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPDestroy(&mglevels[i]->smoothd));
5109566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&mglevels[i]->smoothu));
5119566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&mglevels[i]->cr));
5129566063dSJacob Faibussowitsch       PetscCall(PetscFree(mglevels[i]));
513a06653b4SBarry Smith     }
5149566063dSJacob Faibussowitsch     PetscCall(PetscFree(mg->levels));
515a06653b4SBarry Smith   }
5169566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
5179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL));
5189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL));
5192b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", NULL));
5202b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", NULL));
5212b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", NULL));
5222b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL));
5232b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL));
5242b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", NULL));
5252b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", NULL));
5262b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", NULL));
5272b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", NULL));
5282b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", NULL));
5292b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", NULL));
5303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
531f3fbd535SBarry Smith }
532f3fbd535SBarry Smith 
533f3fbd535SBarry Smith /*
534f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
535f3fbd535SBarry Smith              or full cycle of multigrid.
536f3fbd535SBarry Smith 
537f3fbd535SBarry Smith   Note:
538f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
539f3fbd535SBarry Smith */
540d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_MG_Internal(PC pc, Vec b, Vec x, Mat B, Mat X, PetscBool transpose)
541d71ae5a4SJacob Faibussowitsch {
542f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
543f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
544391689abSStefano Zampini   PC             tpc;
545f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels, i;
54630b0564aSStefano Zampini   PetscBool      changeu, changed, matapp;
547f3fbd535SBarry Smith 
548f3fbd535SBarry Smith   PetscFunctionBegin;
54930b0564aSStefano Zampini   matapp = (PetscBool)(B && X);
5509566063dSJacob Faibussowitsch   if (mg->stageApply) PetscCall(PetscLogStagePush(mg->stageApply));
551e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
552298cc208SBarry Smith   for (i = 0; i < levels; i++) {
553e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
5549566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL));
5559566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A));
556e1d8e5deSBarry Smith     }
557e1d8e5deSBarry Smith   }
558e1d8e5deSBarry Smith 
5599566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc));
5609566063dSJacob Faibussowitsch   PetscCall(PCPreSolveChangeRHS(tpc, &changed));
5619566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc));
5629566063dSJacob Faibussowitsch   PetscCall(PCPreSolveChangeRHS(tpc, &changeu));
563391689abSStefano Zampini   if (!changeu && !changed) {
56430b0564aSStefano Zampini     if (matapp) {
5659566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[levels - 1]->B));
56630b0564aSStefano Zampini       mglevels[levels - 1]->B = B;
56730b0564aSStefano Zampini     } else {
5689566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[levels - 1]->b));
569f3fbd535SBarry Smith       mglevels[levels - 1]->b = b;
57030b0564aSStefano Zampini     }
571391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
57230b0564aSStefano Zampini     if (matapp) {
57330b0564aSStefano Zampini       if (mglevels[levels - 1]->B) {
57430b0564aSStefano Zampini         PetscInt  N1, N2;
57530b0564aSStefano Zampini         PetscBool flg;
57630b0564aSStefano Zampini 
5779566063dSJacob Faibussowitsch         PetscCall(MatGetSize(mglevels[levels - 1]->B, NULL, &N1));
5789566063dSJacob Faibussowitsch         PetscCall(MatGetSize(B, NULL, &N2));
5799566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels - 1]->B, ((PetscObject)B)->type_name, &flg));
58048a46eb9SPierre Jolivet         if (N1 != N2 || !flg) PetscCall(MatDestroy(&mglevels[levels - 1]->B));
58130b0564aSStefano Zampini       }
58230b0564aSStefano Zampini       if (!mglevels[levels - 1]->B) {
5839566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &mglevels[levels - 1]->B));
58430b0564aSStefano Zampini       } else {
5859566063dSJacob Faibussowitsch         PetscCall(MatCopy(B, mglevels[levels - 1]->B, SAME_NONZERO_PATTERN));
58630b0564aSStefano Zampini       }
58730b0564aSStefano Zampini     } else {
588391689abSStefano Zampini       if (!mglevels[levels - 1]->b) {
589391689abSStefano Zampini         Vec *vec;
590391689abSStefano Zampini 
5919566063dSJacob Faibussowitsch         PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL));
592391689abSStefano Zampini         mglevels[levels - 1]->b = *vec;
5939566063dSJacob Faibussowitsch         PetscCall(PetscFree(vec));
594391689abSStefano Zampini       }
5959566063dSJacob Faibussowitsch       PetscCall(VecCopy(b, mglevels[levels - 1]->b));
596391689abSStefano Zampini     }
59730b0564aSStefano Zampini   }
5989371c9d4SSatish Balay   if (matapp) {
5999371c9d4SSatish Balay     mglevels[levels - 1]->X = X;
6009371c9d4SSatish Balay   } else {
6019371c9d4SSatish Balay     mglevels[levels - 1]->x = x;
6029371c9d4SSatish Balay   }
60330b0564aSStefano Zampini 
60430b0564aSStefano Zampini   /* If coarser Xs are present, it means we have already block applied the PC at least once
60530b0564aSStefano Zampini      Reset operators if sizes/type do no match */
60630b0564aSStefano Zampini   if (matapp && levels > 1 && mglevels[levels - 2]->X) {
60730b0564aSStefano Zampini     PetscInt  Xc, Bc;
60830b0564aSStefano Zampini     PetscBool flg;
60930b0564aSStefano Zampini 
6109566063dSJacob Faibussowitsch     PetscCall(MatGetSize(mglevels[levels - 2]->X, NULL, &Xc));
6119566063dSJacob Faibussowitsch     PetscCall(MatGetSize(mglevels[levels - 1]->B, NULL, &Bc));
6129566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels - 2]->X, ((PetscObject)mglevels[levels - 1]->X)->type_name, &flg));
61330b0564aSStefano Zampini     if (Xc != Bc || !flg) {
6149566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[levels - 1]->R));
61530b0564aSStefano Zampini       for (i = 0; i < levels - 1; i++) {
6169566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&mglevels[i]->R));
6179566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&mglevels[i]->B));
6189566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&mglevels[i]->X));
61930b0564aSStefano Zampini       }
62030b0564aSStefano Zampini     }
62130b0564aSStefano Zampini   }
622391689abSStefano Zampini 
62331567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
6249566063dSJacob Faibussowitsch     if (matapp) PetscCall(MatZeroEntries(X));
6259566063dSJacob Faibussowitsch     else PetscCall(VecZeroEntries(x));
62648a46eb9SPierre Jolivet     for (i = 0; i < mg->cyclesperpcapply; i++) PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, transpose, matapp, NULL));
6272fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
6289566063dSJacob Faibussowitsch     PetscCall(PCMGACycle_Private(pc, mglevels, transpose, matapp));
6292fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
6309566063dSJacob Faibussowitsch     PetscCall(PCMGKCycle_Private(pc, mglevels, transpose, matapp));
6312fa5cd67SKarl Rupp   } else {
6329566063dSJacob Faibussowitsch     PetscCall(PCMGFCycle_Private(pc, mglevels, transpose, matapp));
633f3fbd535SBarry Smith   }
6349566063dSJacob Faibussowitsch   if (mg->stageApply) PetscCall(PetscLogStagePop());
63530b0564aSStefano Zampini   if (!changeu && !changed) {
6369371c9d4SSatish Balay     if (matapp) {
6379371c9d4SSatish Balay       mglevels[levels - 1]->B = NULL;
6389371c9d4SSatish Balay     } else {
6399371c9d4SSatish Balay       mglevels[levels - 1]->b = NULL;
6409371c9d4SSatish Balay     }
64130b0564aSStefano Zampini   }
6423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
643f3fbd535SBarry Smith }
644f3fbd535SBarry Smith 
645d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_MG(PC pc, Vec b, Vec x)
646d71ae5a4SJacob Faibussowitsch {
647fcb023d4SJed Brown   PetscFunctionBegin;
6489566063dSJacob Faibussowitsch   PetscCall(PCApply_MG_Internal(pc, b, x, NULL, NULL, PETSC_FALSE));
6493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
650fcb023d4SJed Brown }
651fcb023d4SJed Brown 
652d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_MG(PC pc, Vec b, Vec x)
653d71ae5a4SJacob Faibussowitsch {
654fcb023d4SJed Brown   PetscFunctionBegin;
6559566063dSJacob Faibussowitsch   PetscCall(PCApply_MG_Internal(pc, b, x, NULL, NULL, PETSC_TRUE));
6563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
65730b0564aSStefano Zampini }
65830b0564aSStefano Zampini 
659d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_MG(PC pc, Mat b, Mat x)
660d71ae5a4SJacob Faibussowitsch {
66130b0564aSStefano Zampini   PetscFunctionBegin;
6629566063dSJacob Faibussowitsch   PetscCall(PCApply_MG_Internal(pc, NULL, NULL, b, x, PETSC_FALSE));
6633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
664fcb023d4SJed Brown }
665f3fbd535SBarry Smith 
666d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFromOptions_MG(PC pc, PetscOptionItems *PetscOptionsObject)
667d71ae5a4SJacob Faibussowitsch {
668710315b6SLawrence Mitchell   PetscInt            levels, cycles;
669f3b08a26SMatthew G. Knepley   PetscBool           flg, flg2;
670f3fbd535SBarry Smith   PC_MG              *mg = (PC_MG *)pc->data;
6713d3eaba7SBarry Smith   PC_MG_Levels      **mglevels;
672f3fbd535SBarry Smith   PCMGType            mgtype;
673f3fbd535SBarry Smith   PCMGCycleType       mgctype;
6742134b1e4SBarry Smith   PCMGGalerkinType    gtype;
6752b3cbbdaSStefano Zampini   PCMGCoarseSpaceType coarseSpaceType;
676f3fbd535SBarry Smith 
677f3fbd535SBarry Smith   PetscFunctionBegin;
6781d743356SStefano Zampini   levels = PetscMax(mg->nlevels, 1);
679d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Multigrid options");
6809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_mg_levels", "Number of Levels", "PCMGSetLevels", levels, &levels, &flg));
6811a97d7b6SStefano Zampini   if (!flg && !mg->levels && pc->dm) {
6829566063dSJacob Faibussowitsch     PetscCall(DMGetRefineLevel(pc->dm, &levels));
683eb3f98d2SBarry Smith     levels++;
6843aeaf226SBarry Smith     mg->usedmfornumberoflevels = PETSC_TRUE;
685eb3f98d2SBarry Smith   }
6869566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc, levels, NULL));
6873aeaf226SBarry Smith   mglevels = mg->levels;
6883aeaf226SBarry Smith 
689f3fbd535SBarry Smith   mgctype = (PCMGCycleType)mglevels[0]->cycles;
6909566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_mg_cycle_type", "V cycle or for W-cycle", "PCMGSetCycleType", PCMGCycleTypes, (PetscEnum)mgctype, (PetscEnum *)&mgctype, &flg));
6911baa6e33SBarry Smith   if (flg) PetscCall(PCMGSetCycleType(pc, mgctype));
6922134b1e4SBarry Smith   gtype = mg->galerkin;
6939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_mg_galerkin", "Use Galerkin process to compute coarser operators", "PCMGSetGalerkin", PCMGGalerkinTypes, (PetscEnum)gtype, (PetscEnum *)&gtype, &flg));
6941baa6e33SBarry Smith   if (flg) PetscCall(PCMGSetGalerkin(pc, gtype));
6952b3cbbdaSStefano Zampini   coarseSpaceType = mg->coarseSpaceType;
6962b3cbbdaSStefano 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));
6972b3cbbdaSStefano Zampini   if (flg) PetscCall(PCMGSetAdaptCoarseSpaceType(pc, coarseSpaceType));
6989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_mg_adapt_interp_n", "Size of the coarse space for adaptive interpolation", "PCMGSetCoarseSpace", mg->Nc, &mg->Nc, &flg));
6999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_mg_mesp_monitor", "Monitor the multilevel eigensolver", "PCMGSetAdaptInterpolation", PETSC_FALSE, &mg->mespMonitor, &flg));
70041b6fd38SMatthew G. Knepley   flg2 = PETSC_FALSE;
7019566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_mg_adapt_cr", "Monitor coarse space quality using Compatible Relaxation (CR)", "PCMGSetAdaptCR", PETSC_FALSE, &flg2, &flg));
7029566063dSJacob Faibussowitsch   if (flg) PetscCall(PCMGSetAdaptCR(pc, flg2));
703f56b1265SBarry Smith   flg = PETSC_FALSE;
7049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_mg_distinct_smoothup", "Create separate smoothup KSP and append the prefix _up", "PCMGSetDistinctSmoothUp", PETSC_FALSE, &flg, NULL));
7051baa6e33SBarry Smith   if (flg) PetscCall(PCMGSetDistinctSmoothUp(pc));
70631567311SBarry Smith   mgtype = mg->am;
7079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_mg_type", "Multigrid type", "PCMGSetType", PCMGTypes, (PetscEnum)mgtype, (PetscEnum *)&mgtype, &flg));
7081baa6e33SBarry Smith   if (flg) PetscCall(PCMGSetType(pc, mgtype));
70931567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
7109566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_mg_multiplicative_cycles", "Number of cycles for each preconditioner step", "PCMGMultiplicativeSetCycles", mg->cyclesperpcapply, &cycles, &flg));
7111baa6e33SBarry Smith     if (flg) PetscCall(PCMGMultiplicativeSetCycles(pc, cycles));
712f3fbd535SBarry Smith   }
713f3fbd535SBarry Smith   flg = PETSC_FALSE;
7149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_mg_log", "Log times for each multigrid level", "None", flg, &flg, NULL));
715f3fbd535SBarry Smith   if (flg) {
716f3fbd535SBarry Smith     PetscInt i;
717f3fbd535SBarry Smith     char     eventname[128];
7181a97d7b6SStefano Zampini 
719f3fbd535SBarry Smith     levels = mglevels[0]->levels;
720f3fbd535SBarry Smith     for (i = 0; i < levels; i++) {
721a364092eSJacob Faibussowitsch       PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGSetup Level %d", (int)i));
7229566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventsmoothsetup));
723a364092eSJacob Faibussowitsch       PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGSmooth Level %d", (int)i));
7249566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventsmoothsolve));
725f3fbd535SBarry Smith       if (i) {
726a364092eSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGResid Level %d", (int)i));
7279566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventresidual));
728a364092eSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGInterp Level %d", (int)i));
7299566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventinterprestrict));
730f3fbd535SBarry Smith       }
731f3fbd535SBarry Smith     }
732eec35531SMatthew G Knepley 
7332611ad71SToby Isaac     if (PetscDefined(USE_LOG)) {
7342611ad71SToby Isaac       const char sname[] = "MG Apply";
735eec35531SMatthew G Knepley 
7362611ad71SToby Isaac       PetscCall(PetscLogStageGetId(sname, &mg->stageApply));
7372611ad71SToby Isaac       if (mg->stageApply < 0) PetscCall(PetscLogStageRegister(sname, &mg->stageApply));
738eec35531SMatthew G Knepley     }
739f3fbd535SBarry Smith   }
740d0609cedSBarry Smith   PetscOptionsHeadEnd();
741f3b08a26SMatthew G. Knepley   /* Check option consistency */
7429566063dSJacob Faibussowitsch   PetscCall(PCMGGetGalerkin(pc, &gtype));
7439566063dSJacob Faibussowitsch   PetscCall(PCMGGetAdaptInterpolation(pc, &flg));
74408401ef6SPierre Jolivet   PetscCheck(!flg || !(gtype >= PC_MG_GALERKIN_NONE), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
7453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
746f3fbd535SBarry Smith }
747f3fbd535SBarry Smith 
7480a545947SLisandro Dalcin const char *const PCMGTypes[]            = {"MULTIPLICATIVE", "ADDITIVE", "FULL", "KASKADE", "PCMGType", "PC_MG", NULL};
7490a545947SLisandro Dalcin const char *const PCMGCycleTypes[]       = {"invalid", "v", "w", "PCMGCycleType", "PC_MG_CYCLE", NULL};
7500a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[]    = {"both", "pmat", "mat", "none", "external", "PCMGGalerkinType", "PC_MG_GALERKIN", NULL};
7512b3cbbdaSStefano Zampini const char *const PCMGCoarseSpaceTypes[] = {"none", "polynomial", "harmonic", "eigenvector", "generalized_eigenvector", "gdsw", "PCMGCoarseSpaceType", "PCMG_ADAPT_NONE", NULL};
752f3fbd535SBarry Smith 
7539804daf3SBarry Smith #include <petscdraw.h>
754d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView_MG(PC pc, PetscViewer viewer)
755d71ae5a4SJacob Faibussowitsch {
756f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
757f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
758e3deeeafSJed Brown   PetscInt       levels   = mglevels ? mglevels[0]->levels : 0, i;
759219b1045SBarry Smith   PetscBool      iascii, isbinary, isdraw;
760f3fbd535SBarry Smith 
761f3fbd535SBarry Smith   PetscFunctionBegin;
7629566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
7639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
7649566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
765f3fbd535SBarry Smith   if (iascii) {
766e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
76763a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  type is %s, levels=%" PetscInt_FMT " cycles=%s\n", PCMGTypes[mg->am], levels, cyclename));
76848a46eb9SPierre Jolivet     if (mg->am == PC_MG_MULTIPLICATIVE) PetscCall(PetscViewerASCIIPrintf(viewer, "    Cycles per PCApply=%" PetscInt_FMT "\n", mg->cyclesperpcapply));
7692134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
7709566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using Galerkin computed coarse grid matrices\n"));
7712134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
7729566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using Galerkin computed coarse grid matrices for pmat\n"));
7732134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
7749566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using Galerkin computed coarse grid matrices for mat\n"));
7752134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
7769566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using externally compute Galerkin coarse grid matrices\n"));
7774f66f45eSBarry Smith     } else {
7789566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "    Not using Galerkin computed coarse grid matrices\n"));
779f3fbd535SBarry Smith     }
7801baa6e33SBarry Smith     if (mg->view) PetscCall((*mg->view)(pc, viewer));
781f3fbd535SBarry Smith     for (i = 0; i < levels; i++) {
78263a3b9bcSJacob Faibussowitsch       if (i) {
78363a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n", i));
784f3fbd535SBarry Smith       } else {
78563a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n", i));
786f3fbd535SBarry Smith       }
7879566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
7889566063dSJacob Faibussowitsch       PetscCall(KSPView(mglevels[i]->smoothd, viewer));
7899566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
790f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
7919566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) same as down solver (pre-smoother)\n"));
792f3fbd535SBarry Smith       } else if (i) {
79363a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n", i));
7949566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
7959566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothu, viewer));
7969566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
797f3fbd535SBarry Smith       }
79841b6fd38SMatthew G. Knepley       if (i && mglevels[i]->cr) {
79963a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "CR solver on level %" PetscInt_FMT " -------------------------------\n", i));
8009566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
8019566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->cr, viewer));
8029566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
80341b6fd38SMatthew G. Knepley       }
804f3fbd535SBarry Smith     }
8055b0b0462SBarry Smith   } else if (isbinary) {
8065b0b0462SBarry Smith     for (i = levels - 1; i >= 0; i--) {
8079566063dSJacob Faibussowitsch       PetscCall(KSPView(mglevels[i]->smoothd, viewer));
80848a46eb9SPierre Jolivet       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPView(mglevels[i]->smoothu, viewer));
8095b0b0462SBarry Smith     }
810219b1045SBarry Smith   } else if (isdraw) {
811219b1045SBarry Smith     PetscDraw draw;
812b4375e8dSPeter Brune     PetscReal x, w, y, bottom, th;
8139566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
8149566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
8159566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringGetSize(draw, NULL, &th));
816219b1045SBarry Smith     bottom = y - th;
817219b1045SBarry Smith     for (i = levels - 1; i >= 0; i--) {
818b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
8199566063dSJacob Faibussowitsch         PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
8209566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothd, viewer));
8219566063dSJacob Faibussowitsch         PetscCall(PetscDrawPopCurrentPoint(draw));
822b4375e8dSPeter Brune       } else {
823b4375e8dSPeter Brune         w = 0.5 * PetscMin(1.0 - x, x);
8249566063dSJacob Faibussowitsch         PetscCall(PetscDrawPushCurrentPoint(draw, x + w, bottom));
8259566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothd, viewer));
8269566063dSJacob Faibussowitsch         PetscCall(PetscDrawPopCurrentPoint(draw));
8279566063dSJacob Faibussowitsch         PetscCall(PetscDrawPushCurrentPoint(draw, x - w, bottom));
8289566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothu, viewer));
8299566063dSJacob Faibussowitsch         PetscCall(PetscDrawPopCurrentPoint(draw));
830b4375e8dSPeter Brune       }
8319566063dSJacob Faibussowitsch       PetscCall(PetscDrawGetBoundingBox(draw, NULL, &bottom, NULL, NULL));
8321cd381d6SBarry Smith       bottom -= th;
833219b1045SBarry Smith     }
8345b0b0462SBarry Smith   }
8353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
836f3fbd535SBarry Smith }
837f3fbd535SBarry Smith 
838af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
839cab2e9ccSBarry Smith 
840f3fbd535SBarry Smith /*
841f3fbd535SBarry Smith     Calls setup for the KSP on each level
842f3fbd535SBarry Smith */
843d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp_MG(PC pc)
844d71ae5a4SJacob Faibussowitsch {
845f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
846f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
8471a97d7b6SStefano Zampini   PetscInt       i, n;
84898e04984SBarry Smith   PC             cpc;
8493db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE, opsset, use_amat, missinginterpolate = PETSC_FALSE;
850f3fbd535SBarry Smith   Mat            dA, dB;
851f3fbd535SBarry Smith   Vec            tvec;
852218a07d4SBarry Smith   DM            *dms;
8530a545947SLisandro Dalcin   PetscViewer    viewer = NULL;
85441b6fd38SMatthew G. Knepley   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE;
8552b3cbbdaSStefano Zampini   PetscBool      adaptInterpolation = mg->adaptInterpolation;
856f3fbd535SBarry Smith 
857f3fbd535SBarry Smith   PetscFunctionBegin;
85828b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels with PCMGSetLevels() before setting up");
8591a97d7b6SStefano Zampini   n = mglevels[0]->levels;
86001bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
8613aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
8623aeaf226SBarry Smith     PetscInt levels;
8639566063dSJacob Faibussowitsch     PetscCall(DMGetRefineLevel(pc->dm, &levels));
8643aeaf226SBarry Smith     levels++;
8653aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
8669566063dSJacob Faibussowitsch       PetscCall(PCMGSetLevels(pc, levels, NULL));
8673aeaf226SBarry Smith       n = levels;
8689566063dSJacob Faibussowitsch       PetscCall(PCSetFromOptions(pc)); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
8693aeaf226SBarry Smith       mglevels = mg->levels;
8703aeaf226SBarry Smith     }
8713aeaf226SBarry Smith   }
8729566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[0]->smoothd, &cpc));
8733aeaf226SBarry Smith 
874f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
875f3fbd535SBarry Smith   /* so use those from global PC */
876f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
8779566063dSJacob Faibussowitsch   PetscCall(KSPGetOperatorsSet(mglevels[n - 1]->smoothd, NULL, &opsset));
87898e04984SBarry Smith   if (opsset) {
87998e04984SBarry Smith     Mat mmat;
8809566063dSJacob Faibussowitsch     PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, NULL, &mmat));
88198e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
88298e04984SBarry Smith   }
8834a5f07a5SMark F. Adams 
88441b6fd38SMatthew G. Knepley   /* Create CR solvers */
8859566063dSJacob Faibussowitsch   PetscCall(PCMGGetAdaptCR(pc, &doCR));
88641b6fd38SMatthew G. Knepley   if (doCR) {
88741b6fd38SMatthew G. Knepley     const char *prefix;
88841b6fd38SMatthew G. Knepley 
8899566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
89041b6fd38SMatthew G. Knepley     for (i = 1; i < n; ++i) {
89141b6fd38SMatthew G. Knepley       PC   ipc, cr;
89241b6fd38SMatthew G. Knepley       char crprefix[128];
89341b6fd38SMatthew G. Knepley 
8949566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &mglevels[i]->cr));
8953821be0aSBarry Smith       PetscCall(KSPSetNestLevel(mglevels[i]->cr, pc->kspnestlevel));
8969566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE));
8979566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->cr, (PetscObject)pc, n - i));
8989566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(mglevels[i]->cr, prefix));
8999566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level));
9009566063dSJacob Faibussowitsch       PetscCall(KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV));
9019566063dSJacob Faibussowitsch       PetscCall(KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL));
9029566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED));
9039566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(mglevels[i]->cr, &ipc));
90441b6fd38SMatthew G. Knepley 
9059566063dSJacob Faibussowitsch       PetscCall(PCSetType(ipc, PCCOMPOSITE));
9069566063dSJacob Faibussowitsch       PetscCall(PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE));
9079566063dSJacob Faibussowitsch       PetscCall(PCCompositeAddPCType(ipc, PCSOR));
9089566063dSJacob Faibussowitsch       PetscCall(CreateCR_Private(pc, i, &cr));
9099566063dSJacob Faibussowitsch       PetscCall(PCCompositeAddPC(ipc, cr));
9109566063dSJacob Faibussowitsch       PetscCall(PCDestroy(&cr));
91141b6fd38SMatthew G. Knepley 
9129566063dSJacob Faibussowitsch       PetscCall(KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd));
9139566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE));
9149566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int)i));
9159566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix));
91641b6fd38SMatthew G. Knepley     }
91741b6fd38SMatthew G. Knepley   }
91841b6fd38SMatthew G. Knepley 
9194a5f07a5SMark F. Adams   if (!opsset) {
9209566063dSJacob Faibussowitsch     PetscCall(PCGetUseAmat(pc, &use_amat));
9214a5f07a5SMark F. Adams     if (use_amat) {
9229566063dSJacob 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"));
9239566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->mat, pc->pmat));
92469858f1bSStefano Zampini     } else {
9259566063dSJacob 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"));
9269566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->pmat, pc->pmat));
9274a5f07a5SMark F. Adams     }
9284a5f07a5SMark F. Adams   }
929f3fbd535SBarry Smith 
9309c7ffb39SBarry Smith   for (i = n - 1; i > 0; i--) {
9319c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
9329c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
9332b3cbbdaSStefano Zampini       break;
9349c7ffb39SBarry Smith     }
9359c7ffb39SBarry Smith   }
9362134b1e4SBarry Smith 
9379566063dSJacob Faibussowitsch   PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, &dA, &dB));
9382134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
9392b3cbbdaSStefano Zampini   if (mg->galerkin == PC_MG_GALERKIN_NONE || ((mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_MAT) && !dAeqdB)) {
9402134b1e4SBarry Smith     needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
9412134b1e4SBarry Smith   }
9422134b1e4SBarry Smith 
9432b3cbbdaSStefano Zampini   if (pc->dm && !pc->setupcalled) {
9442b3cbbdaSStefano Zampini     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
9452b3cbbdaSStefano Zampini     PetscCall(KSPSetDM(mglevels[n - 1]->smoothd, pc->dm));
9462b3cbbdaSStefano Zampini     PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothd, PETSC_FALSE));
9472b3cbbdaSStefano Zampini     if (mglevels[n - 1]->smoothd != mglevels[n - 1]->smoothu) {
9482b3cbbdaSStefano Zampini       PetscCall(KSPSetDM(mglevels[n - 1]->smoothu, pc->dm));
9492b3cbbdaSStefano Zampini       PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothu, PETSC_FALSE));
9502b3cbbdaSStefano Zampini     }
9512b3cbbdaSStefano Zampini     if (mglevels[n - 1]->cr) {
9522b3cbbdaSStefano Zampini       PetscCall(KSPSetDM(mglevels[n - 1]->cr, pc->dm));
9532b3cbbdaSStefano Zampini       PetscCall(KSPSetDMActive(mglevels[n - 1]->cr, PETSC_FALSE));
9542b3cbbdaSStefano Zampini     }
9552b3cbbdaSStefano Zampini   }
9562b3cbbdaSStefano Zampini 
9579c7ffb39SBarry Smith   /*
9589c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
9592b3cbbdaSStefano Zampini    Skipping for externally managed hierarchy (such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
9609c7ffb39SBarry Smith   */
96132fe681dSStefano Zampini   if (missinginterpolate && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
96232fe681dSStefano Zampini     /* first see if we can compute a coarse space */
96332fe681dSStefano Zampini     if (mg->coarseSpaceType == PCMG_ADAPT_GDSW) {
96432fe681dSStefano Zampini       for (i = n - 2; i > -1; i--) {
96532fe681dSStefano Zampini         if (!mglevels[i + 1]->restrct && !mglevels[i + 1]->interpolate) {
96632fe681dSStefano Zampini           PetscCall(PCMGComputeCoarseSpace_Internal(pc, i + 1, mg->coarseSpaceType, mg->Nc, NULL, &mglevels[i + 1]->coarseSpace));
96732fe681dSStefano Zampini           PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->coarseSpace));
96832fe681dSStefano Zampini         }
96932fe681dSStefano Zampini       }
97032fe681dSStefano Zampini     } else { /* construct the interpolation from the DMs */
9712e499ae9SBarry Smith       Mat p;
97273250ac0SBarry Smith       Vec rscale;
9739566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &dms));
974218a07d4SBarry Smith       dms[n - 1] = pc->dm;
975ef1267afSMatthew G. Knepley       /* Separately create them so we do not get DMKSP interference between levels */
9769566063dSJacob Faibussowitsch       for (i = n - 2; i > -1; i--) PetscCall(DMCoarsen(dms[i + 1], MPI_COMM_NULL, &dms[i]));
977218a07d4SBarry Smith       for (i = n - 2; i > -1; i--) {
978942e3340SBarry Smith         DMKSP     kdm;
979eab52d2bSLawrence Mitchell         PetscBool dmhasrestrict, dmhasinject;
9809566063dSJacob Faibussowitsch         PetscCall(KSPSetDM(mglevels[i]->smoothd, dms[i]));
9819566063dSJacob Faibussowitsch         if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothd, PETSC_FALSE));
982c27ee7a3SPatrick Farrell         if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
9839566063dSJacob Faibussowitsch           PetscCall(KSPSetDM(mglevels[i]->smoothu, dms[i]));
9849566063dSJacob Faibussowitsch           if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothu, PETSC_FALSE));
985c27ee7a3SPatrick Farrell         }
98641b6fd38SMatthew G. Knepley         if (mglevels[i]->cr) {
9879566063dSJacob Faibussowitsch           PetscCall(KSPSetDM(mglevels[i]->cr, dms[i]));
9889566063dSJacob Faibussowitsch           if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->cr, PETSC_FALSE));
98941b6fd38SMatthew G. Knepley         }
9909566063dSJacob Faibussowitsch         PetscCall(DMGetDMKSPWrite(dms[i], &kdm));
991d1a3e677SJed Brown         /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
992d1a3e677SJed Brown          * a bitwise OR of computing the matrix, RHS, and initial iterate. */
9930298fd71SBarry Smith         kdm->ops->computerhs = NULL;
9940298fd71SBarry Smith         kdm->rhsctx          = NULL;
99524c3aa18SBarry Smith         if (!mglevels[i + 1]->interpolate) {
9969566063dSJacob Faibussowitsch           PetscCall(DMCreateInterpolation(dms[i], dms[i + 1], &p, &rscale));
9979566063dSJacob Faibussowitsch           PetscCall(PCMGSetInterpolation(pc, i + 1, p));
9989566063dSJacob Faibussowitsch           if (rscale) PetscCall(PCMGSetRScale(pc, i + 1, rscale));
9999566063dSJacob Faibussowitsch           PetscCall(VecDestroy(&rscale));
10009566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&p));
1001218a07d4SBarry Smith         }
10029566063dSJacob Faibussowitsch         PetscCall(DMHasCreateRestriction(dms[i], &dmhasrestrict));
10033ad4599aSBarry Smith         if (dmhasrestrict && !mglevels[i + 1]->restrct) {
10049566063dSJacob Faibussowitsch           PetscCall(DMCreateRestriction(dms[i], dms[i + 1], &p));
10059566063dSJacob Faibussowitsch           PetscCall(PCMGSetRestriction(pc, i + 1, p));
10069566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&p));
10073ad4599aSBarry Smith         }
10089566063dSJacob Faibussowitsch         PetscCall(DMHasCreateInjection(dms[i], &dmhasinject));
1009eab52d2bSLawrence Mitchell         if (dmhasinject && !mglevels[i + 1]->inject) {
10109566063dSJacob Faibussowitsch           PetscCall(DMCreateInjection(dms[i], dms[i + 1], &p));
10119566063dSJacob Faibussowitsch           PetscCall(PCMGSetInjection(pc, i + 1, p));
10129566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&p));
1013eab52d2bSLawrence Mitchell         }
101424c3aa18SBarry Smith       }
10152d2b81a6SBarry Smith 
10169566063dSJacob Faibussowitsch       for (i = n - 2; i > -1; i--) PetscCall(DMDestroy(&dms[i]));
10179566063dSJacob Faibussowitsch       PetscCall(PetscFree(dms));
1018ce4cda84SJed Brown     }
101932fe681dSStefano Zampini   }
1020cab2e9ccSBarry Smith 
10212134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
10222134b1e4SBarry Smith     Mat       A, B;
10232134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE, doB = PETSC_FALSE;
10242134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
10252134b1e4SBarry Smith 
10262b3cbbdaSStefano Zampini     if (mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_BOTH) doB = PETSC_TRUE;
10272b3cbbdaSStefano Zampini     if (mg->galerkin == PC_MG_GALERKIN_MAT || (mg->galerkin == PC_MG_GALERKIN_BOTH && dA != dB)) doA = PETSC_TRUE;
10282134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
1029f3fbd535SBarry Smith     for (i = n - 2; i > -1; i--) {
10307827d75bSBarry 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");
103148a46eb9SPierre Jolivet       if (!mglevels[i + 1]->interpolate) PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->restrct));
103248a46eb9SPierre Jolivet       if (!mglevels[i + 1]->restrct) PetscCall(PCMGSetRestriction(pc, i + 1, mglevels[i + 1]->interpolate));
103348a46eb9SPierre Jolivet       if (reuse == MAT_REUSE_MATRIX) PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, &B));
103448a46eb9SPierre Jolivet       if (doA) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dA, mglevels[i + 1]->interpolate, reuse, 1.0, &A));
103548a46eb9SPierre Jolivet       if (doB) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dB, mglevels[i + 1]->interpolate, reuse, 1.0, &B));
10362134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
10372134b1e4SBarry Smith       if (!doA && dAeqdB) {
10389566063dSJacob Faibussowitsch         if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)B));
10392134b1e4SBarry Smith         A = B;
10402134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
10419566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, NULL));
10429566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)A));
1043b9367914SBarry Smith       }
10442134b1e4SBarry Smith       if (!doB && dAeqdB) {
10459566063dSJacob Faibussowitsch         if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)A));
10462134b1e4SBarry Smith         B = A;
10472134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
10489566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd, NULL, &B));
10499566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)B));
10507d537d92SJed Brown       }
10512134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
10529566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[i]->smoothd, A, B));
10539566063dSJacob Faibussowitsch         PetscCall(PetscObjectDereference((PetscObject)A));
10549566063dSJacob Faibussowitsch         PetscCall(PetscObjectDereference((PetscObject)B));
10552134b1e4SBarry Smith       }
10562134b1e4SBarry Smith       dA = A;
1057f3fbd535SBarry Smith       dB = B;
1058f3fbd535SBarry Smith     }
1059f3fbd535SBarry Smith   }
1060f3b08a26SMatthew G. Knepley 
1061f3b08a26SMatthew G. Knepley   /* Adapt interpolation matrices */
10622b3cbbdaSStefano Zampini   if (adaptInterpolation) {
1063f3b08a26SMatthew G. Knepley     for (i = 0; i < n; ++i) {
106448a46eb9SPierre Jolivet       if (!mglevels[i]->coarseSpace) PetscCall(PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i - 1]->coarseSpace, &mglevels[i]->coarseSpace));
10652b3cbbdaSStefano Zampini       if (i) PetscCall(PCMGAdaptInterpolator_Internal(pc, i, mglevels[i - 1]->smoothu, mglevels[i]->smoothu, mglevels[i - 1]->coarseSpace, mglevels[i]->coarseSpace));
1066f3b08a26SMatthew G. Knepley     }
106748a46eb9SPierre Jolivet     for (i = n - 2; i > -1; --i) PetscCall(PCMGRecomputeLevelOperators_Internal(pc, i));
1068f3b08a26SMatthew G. Knepley   }
1069f3b08a26SMatthew G. Knepley 
10702134b1e4SBarry Smith   if (needRestricts && pc->dm) {
1071caa4e7f2SJed Brown     for (i = n - 2; i >= 0; i--) {
1072ccceb50cSJed Brown       DM  dmfine, dmcoarse;
1073caa4e7f2SJed Brown       Mat Restrict, Inject;
1074caa4e7f2SJed Brown       Vec rscale;
10759566063dSJacob Faibussowitsch       PetscCall(KSPGetDM(mglevels[i + 1]->smoothd, &dmfine));
10769566063dSJacob Faibussowitsch       PetscCall(KSPGetDM(mglevels[i]->smoothd, &dmcoarse));
10779566063dSJacob Faibussowitsch       PetscCall(PCMGGetRestriction(pc, i + 1, &Restrict));
10789566063dSJacob Faibussowitsch       PetscCall(PCMGGetRScale(pc, i + 1, &rscale));
10799566063dSJacob Faibussowitsch       PetscCall(PCMGGetInjection(pc, i + 1, &Inject));
10809566063dSJacob Faibussowitsch       PetscCall(DMRestrict(dmfine, Restrict, rscale, Inject, dmcoarse));
1081caa4e7f2SJed Brown     }
1082caa4e7f2SJed Brown   }
1083f3fbd535SBarry Smith 
1084f3fbd535SBarry Smith   if (!pc->setupcalled) {
108548a46eb9SPierre Jolivet     for (i = 0; i < n; i++) PetscCall(KSPSetFromOptions(mglevels[i]->smoothd));
1086f3fbd535SBarry Smith     for (i = 1; i < n; i++) {
108748a46eb9SPierre Jolivet       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) PetscCall(KSPSetFromOptions(mglevels[i]->smoothu));
108848a46eb9SPierre Jolivet       if (mglevels[i]->cr) PetscCall(KSPSetFromOptions(mglevels[i]->cr));
1089f3fbd535SBarry Smith     }
109015229ffcSPierre Jolivet     /* insure that if either interpolation or restriction is set the other one is set */
1091f3fbd535SBarry Smith     for (i = 1; i < n; i++) {
10929566063dSJacob Faibussowitsch       PetscCall(PCMGGetInterpolation(pc, i, NULL));
10939566063dSJacob Faibussowitsch       PetscCall(PCMGGetRestriction(pc, i, NULL));
1094f3fbd535SBarry Smith     }
1095f3fbd535SBarry Smith     for (i = 0; i < n - 1; i++) {
1096f3fbd535SBarry Smith       if (!mglevels[i]->b) {
1097f3fbd535SBarry Smith         Vec *vec;
10989566063dSJacob Faibussowitsch         PetscCall(KSPCreateVecs(mglevels[i]->smoothd, 1, &vec, 0, NULL));
10999566063dSJacob Faibussowitsch         PetscCall(PCMGSetRhs(pc, i, *vec));
11009566063dSJacob Faibussowitsch         PetscCall(VecDestroy(vec));
11019566063dSJacob Faibussowitsch         PetscCall(PetscFree(vec));
1102f3fbd535SBarry Smith       }
1103f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
11049566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b, &tvec));
11059566063dSJacob Faibussowitsch         PetscCall(PCMGSetR(pc, i, tvec));
11069566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&tvec));
1107f3fbd535SBarry Smith       }
1108f3fbd535SBarry Smith       if (!mglevels[i]->x) {
11099566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b, &tvec));
11109566063dSJacob Faibussowitsch         PetscCall(PCMGSetX(pc, i, tvec));
11119566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&tvec));
1112f3fbd535SBarry Smith       }
111341b6fd38SMatthew G. Knepley       if (doCR) {
11149566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crx));
11159566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crb));
11169566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(mglevels[i]->crb));
111741b6fd38SMatthew G. Knepley       }
1118f3fbd535SBarry Smith     }
1119f3fbd535SBarry Smith     if (n != 1 && !mglevels[n - 1]->r) {
1120f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
1121f3fbd535SBarry Smith       Vec *vec;
11229566063dSJacob Faibussowitsch       PetscCall(KSPCreateVecs(mglevels[n - 1]->smoothd, 1, &vec, 0, NULL));
11239566063dSJacob Faibussowitsch       PetscCall(PCMGSetR(pc, n - 1, *vec));
11249566063dSJacob Faibussowitsch       PetscCall(VecDestroy(vec));
11259566063dSJacob Faibussowitsch       PetscCall(PetscFree(vec));
1126f3fbd535SBarry Smith     }
112741b6fd38SMatthew G. Knepley     if (doCR) {
11289566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crx));
11299566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crb));
11309566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(mglevels[n - 1]->crb));
113141b6fd38SMatthew G. Knepley     }
1132f3fbd535SBarry Smith   }
1133f3fbd535SBarry Smith 
113498e04984SBarry Smith   if (pc->dm) {
113598e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
113698e04984SBarry Smith     for (i = 0; i < n - 1; i++) {
113798e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
113898e04984SBarry Smith     }
113998e04984SBarry Smith   }
114008549ed4SJed Brown   // We got here (PCSetUp_MG) because the matrix has changed, which means the smoother needs to be set up again (e.g.,
114108549ed4SJed Brown   // new diagonal for Jacobi). Setting it here allows it to be logged under PCSetUp rather than deep inside a PCApply.
114208549ed4SJed Brown   if (mglevels[n - 1]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[n - 1]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
1143f3fbd535SBarry Smith 
1144f3fbd535SBarry Smith   for (i = 1; i < n; i++) {
11452a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1) {
1146f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
11479566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE));
1148f3fbd535SBarry Smith     }
11499566063dSJacob Faibussowitsch     if (mglevels[i]->cr) PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE));
11509566063dSJacob Faibussowitsch     if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
11519566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(mglevels[i]->smoothd));
1152d4645a75SStefano Zampini     if (mglevels[i]->smoothd->reason) pc->failedreason = PC_SUBPC_ERROR;
11539566063dSJacob Faibussowitsch     if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
1154d42688cbSBarry Smith     if (!mglevels[i]->residual) {
1155d42688cbSBarry Smith       Mat mat;
11569566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL));
11579566063dSJacob Faibussowitsch       PetscCall(PCMGSetResidual(pc, i, PCMGResidualDefault, mat));
1158d42688cbSBarry Smith     }
1159fcb023d4SJed Brown     if (!mglevels[i]->residualtranspose) {
1160fcb023d4SJed Brown       Mat mat;
11619566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL));
11629566063dSJacob Faibussowitsch       PetscCall(PCMGSetResidualTranspose(pc, i, PCMGResidualTransposeDefault, mat));
1163fcb023d4SJed Brown     }
1164f3fbd535SBarry Smith   }
1165f3fbd535SBarry Smith   for (i = 1; i < n; i++) {
1166f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
1167f3fbd535SBarry Smith       Mat downmat, downpmat;
1168f3fbd535SBarry Smith 
1169f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
11709566063dSJacob Faibussowitsch       PetscCall(KSPGetOperatorsSet(mglevels[i]->smoothu, &opsset, NULL));
1171f3fbd535SBarry Smith       if (!opsset) {
11729566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat));
11739566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[i]->smoothu, downmat, downpmat));
1174f3fbd535SBarry Smith       }
1175f3fbd535SBarry Smith 
11769566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothu, PETSC_TRUE));
11779566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
11789566063dSJacob Faibussowitsch       PetscCall(KSPSetUp(mglevels[i]->smoothu));
1179d4645a75SStefano Zampini       if (mglevels[i]->smoothu->reason) pc->failedreason = PC_SUBPC_ERROR;
11809566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
1181f3fbd535SBarry Smith     }
118241b6fd38SMatthew G. Knepley     if (mglevels[i]->cr) {
118341b6fd38SMatthew G. Knepley       Mat downmat, downpmat;
118441b6fd38SMatthew G. Knepley 
118541b6fd38SMatthew G. Knepley       /* check if operators have been set for up, if not use down operators to set them */
11869566063dSJacob Faibussowitsch       PetscCall(KSPGetOperatorsSet(mglevels[i]->cr, &opsset, NULL));
118741b6fd38SMatthew G. Knepley       if (!opsset) {
11889566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat));
11899566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[i]->cr, downmat, downpmat));
119041b6fd38SMatthew G. Knepley       }
119141b6fd38SMatthew G. Knepley 
11929566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE));
11939566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
11949566063dSJacob Faibussowitsch       PetscCall(KSPSetUp(mglevels[i]->cr));
1195d4645a75SStefano Zampini       if (mglevels[i]->cr->reason) pc->failedreason = PC_SUBPC_ERROR;
11969566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
119741b6fd38SMatthew G. Knepley     }
1198f3fbd535SBarry Smith   }
1199f3fbd535SBarry Smith 
12009566063dSJacob Faibussowitsch   if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0));
12019566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(mglevels[0]->smoothd));
1202d4645a75SStefano Zampini   if (mglevels[0]->smoothd->reason) pc->failedreason = PC_SUBPC_ERROR;
12039566063dSJacob Faibussowitsch   if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0));
1204f3fbd535SBarry Smith 
1205f3fbd535SBarry Smith     /*
1206f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
1207e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
1208f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
1209f3fbd535SBarry Smith 
1210f3fbd535SBarry Smith    Only support one or the other at the same time.
1211f3fbd535SBarry Smith   */
1212f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
12139566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_matlab", &dump, NULL));
1214ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
1215f3fbd535SBarry Smith   dump = PETSC_FALSE;
1216f3fbd535SBarry Smith #endif
12179566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_binary", &dump, NULL));
1218ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
1219f3fbd535SBarry Smith 
1220f3fbd535SBarry Smith   if (viewer) {
122148a46eb9SPierre Jolivet     for (i = 1; i < n; i++) PetscCall(MatView(mglevels[i]->restrct, viewer));
1222f3fbd535SBarry Smith     for (i = 0; i < n; i++) {
12239566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(mglevels[i]->smoothd, &pc));
12249566063dSJacob Faibussowitsch       PetscCall(MatView(pc->mat, viewer));
1225f3fbd535SBarry Smith     }
1226f3fbd535SBarry Smith   }
12273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1228f3fbd535SBarry Smith }
1229f3fbd535SBarry Smith 
1230d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
1231d71ae5a4SJacob Faibussowitsch {
123297d33e41SMatthew G. Knepley   PC_MG *mg = (PC_MG *)pc->data;
123397d33e41SMatthew G. Knepley 
123497d33e41SMatthew G. Knepley   PetscFunctionBegin;
123597d33e41SMatthew G. Knepley   *levels = mg->nlevels;
12363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
123797d33e41SMatthew G. Knepley }
123897d33e41SMatthew G. Knepley 
12394b9ad928SBarry Smith /*@
1240f1580f4eSBarry Smith   PCMGGetLevels - Gets the number of levels to use with `PCMG`.
12414b9ad928SBarry Smith 
12424b9ad928SBarry Smith   Not Collective
12434b9ad928SBarry Smith 
12444b9ad928SBarry Smith   Input Parameter:
12454b9ad928SBarry Smith . pc - the preconditioner context
12464b9ad928SBarry Smith 
1247feefa0e1SJacob Faibussowitsch   Output Parameter:
12484b9ad928SBarry Smith . levels - the number of levels
12494b9ad928SBarry Smith 
12504b9ad928SBarry Smith   Level: advanced
12514b9ad928SBarry Smith 
1252562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetLevels()`
12534b9ad928SBarry Smith @*/
1254d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetLevels(PC pc, PetscInt *levels)
1255d71ae5a4SJacob Faibussowitsch {
12564b9ad928SBarry Smith   PetscFunctionBegin;
12570700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12584f572ea9SToby Isaac   PetscAssertPointer(levels, 2);
125997d33e41SMatthew G. Knepley   *levels = 0;
1260cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGGetLevels_C", (PC, PetscInt *), (pc, levels));
12613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12624b9ad928SBarry Smith }
12634b9ad928SBarry Smith 
12644b9ad928SBarry Smith /*@
1265f1580f4eSBarry Smith   PCMGGetGridComplexity - compute operator and grid complexity of the `PCMG` hierarchy
1266e7d4b4cbSMark Adams 
1267e7d4b4cbSMark Adams   Input Parameter:
1268e7d4b4cbSMark Adams . pc - the preconditioner context
1269e7d4b4cbSMark Adams 
127090db8557SMark Adams   Output Parameters:
127190db8557SMark Adams + gc - grid complexity = sum_i(n_i) / n_0
127290db8557SMark Adams - oc - operator complexity = sum_i(nnz_i) / nnz_0
1273e7d4b4cbSMark Adams 
1274e7d4b4cbSMark Adams   Level: advanced
1275e7d4b4cbSMark Adams 
1276f1580f4eSBarry Smith   Note:
1277f1580f4eSBarry Smith   This is often call the operator complexity in multigrid literature
1278f1580f4eSBarry Smith 
1279562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`
1280e7d4b4cbSMark Adams @*/
1281d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc, PetscReal *oc)
1282d71ae5a4SJacob Faibussowitsch {
1283e7d4b4cbSMark Adams   PC_MG         *mg       = (PC_MG *)pc->data;
1284e7d4b4cbSMark Adams   PC_MG_Levels **mglevels = mg->levels;
1285e7d4b4cbSMark Adams   PetscInt       lev, N;
1286e7d4b4cbSMark Adams   PetscLogDouble nnz0 = 0, sgc = 0, soc = 0, n0 = 0;
1287e7d4b4cbSMark Adams   MatInfo        info;
1288e7d4b4cbSMark Adams 
1289e7d4b4cbSMark Adams   PetscFunctionBegin;
129090db8557SMark Adams   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
12914f572ea9SToby Isaac   if (gc) PetscAssertPointer(gc, 2);
12924f572ea9SToby Isaac   if (oc) PetscAssertPointer(oc, 3);
1293e7d4b4cbSMark Adams   if (!pc->setupcalled) {
129490db8557SMark Adams     if (gc) *gc = 0;
129590db8557SMark Adams     if (oc) *oc = 0;
12963ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1297e7d4b4cbSMark Adams   }
1298e7d4b4cbSMark Adams   PetscCheck(mg->nlevels > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MG has no levels");
1299e7d4b4cbSMark Adams   for (lev = 0; lev < mg->nlevels; lev++) {
1300e7d4b4cbSMark Adams     Mat dB;
13019566063dSJacob Faibussowitsch     PetscCall(KSPGetOperators(mglevels[lev]->smoothd, NULL, &dB));
13029566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(dB, MAT_GLOBAL_SUM, &info)); /* global reduction */
13039566063dSJacob Faibussowitsch     PetscCall(MatGetSize(dB, &N, NULL));
1304e7d4b4cbSMark Adams     sgc += N;
1305e7d4b4cbSMark Adams     soc += info.nz_used;
13069371c9d4SSatish Balay     if (lev == mg->nlevels - 1) {
13079371c9d4SSatish Balay       nnz0 = info.nz_used;
13089371c9d4SSatish Balay       n0   = N;
13099371c9d4SSatish Balay     }
1310e7d4b4cbSMark Adams   }
13110fdf79fbSJacob Faibussowitsch   PetscCheck(n0 > 0 && gc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number for grid points on finest level is not available");
13120fdf79fbSJacob Faibussowitsch   *gc = (PetscReal)(sgc / n0);
131390db8557SMark Adams   if (nnz0 > 0 && oc) *oc = (PetscReal)(soc / nnz0);
13143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1315e7d4b4cbSMark Adams }
1316e7d4b4cbSMark Adams 
1317e7d4b4cbSMark Adams /*@
131804c3f3b8SBarry Smith   PCMGSetType - Determines the form of multigrid to use, either
13194b9ad928SBarry Smith   multiplicative, additive, full, or the Kaskade algorithm.
13204b9ad928SBarry Smith 
1321c3339decSBarry Smith   Logically Collective
13224b9ad928SBarry Smith 
13234b9ad928SBarry Smith   Input Parameters:
13244b9ad928SBarry Smith + pc   - the preconditioner context
1325f1580f4eSBarry Smith - form - multigrid form, one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE`
13264b9ad928SBarry Smith 
13274b9ad928SBarry Smith   Options Database Key:
132820f4b53cSBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, additive, full, kaskade
13294b9ad928SBarry Smith 
13304b9ad928SBarry Smith   Level: advanced
13314b9ad928SBarry Smith 
1332562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGGetType()`, `PCMGCycleType`
13334b9ad928SBarry Smith @*/
1334d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetType(PC pc, PCMGType form)
1335d71ae5a4SJacob Faibussowitsch {
1336f3fbd535SBarry Smith   PC_MG *mg = (PC_MG *)pc->data;
13374b9ad928SBarry Smith 
13384b9ad928SBarry Smith   PetscFunctionBegin;
13390700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1340c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc, form, 2);
134131567311SBarry Smith   mg->am = form;
13429dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1343c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
13443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1345c60c7ad4SBarry Smith }
1346c60c7ad4SBarry Smith 
1347c60c7ad4SBarry Smith /*@
1348f1580f4eSBarry Smith   PCMGGetType - Finds the form of multigrid the `PCMG` is using  multiplicative, additive, full, or the Kaskade algorithm.
1349c60c7ad4SBarry Smith 
1350c3339decSBarry Smith   Logically Collective
1351c60c7ad4SBarry Smith 
1352c60c7ad4SBarry Smith   Input Parameter:
1353c60c7ad4SBarry Smith . pc - the preconditioner context
1354c60c7ad4SBarry Smith 
1355c60c7ad4SBarry Smith   Output Parameter:
1356f1580f4eSBarry Smith . type - one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE`, `PCMGCycleType`
1357c60c7ad4SBarry Smith 
1358c60c7ad4SBarry Smith   Level: advanced
1359c60c7ad4SBarry Smith 
1360562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGSetType()`
1361c60c7ad4SBarry Smith @*/
1362d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetType(PC pc, PCMGType *type)
1363d71ae5a4SJacob Faibussowitsch {
1364c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG *)pc->data;
1365c60c7ad4SBarry Smith 
1366c60c7ad4SBarry Smith   PetscFunctionBegin;
1367c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1368c60c7ad4SBarry Smith   *type = mg->am;
13693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13704b9ad928SBarry Smith }
13714b9ad928SBarry Smith 
13724b9ad928SBarry Smith /*@
1373f1580f4eSBarry Smith   PCMGSetCycleType - Sets the type cycles to use.  Use `PCMGSetCycleTypeOnLevel()` for more
13744b9ad928SBarry Smith   complicated cycling.
13754b9ad928SBarry Smith 
1376c3339decSBarry Smith   Logically Collective
13774b9ad928SBarry Smith 
13784b9ad928SBarry Smith   Input Parameters:
1379c2be2410SBarry Smith + pc - the multigrid context
1380f1580f4eSBarry Smith - n  - either `PC_MG_CYCLE_V` or `PC_MG_CYCLE_W`
13814b9ad928SBarry Smith 
13824b9ad928SBarry Smith   Options Database Key:
1383c1cbb1deSBarry Smith . -pc_mg_cycle_type <v,w> - provide the cycle desired
13844b9ad928SBarry Smith 
13854b9ad928SBarry Smith   Level: advanced
13864b9ad928SBarry Smith 
1387562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetCycleTypeOnLevel()`, `PCMGType`, `PCMGCycleType`
13884b9ad928SBarry Smith @*/
1389d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetCycleType(PC pc, PCMGCycleType n)
1390d71ae5a4SJacob Faibussowitsch {
1391f3fbd535SBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
1392f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
139379416396SBarry Smith   PetscInt       i, levels;
13944b9ad928SBarry Smith 
13954b9ad928SBarry Smith   PetscFunctionBegin;
13960700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
139769fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc, n, 2);
139828b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling");
1399f3fbd535SBarry Smith   levels = mglevels[0]->levels;
14002fa5cd67SKarl Rupp   for (i = 0; i < levels; i++) mglevels[i]->cycles = n;
14013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14024b9ad928SBarry Smith }
14034b9ad928SBarry Smith 
14048cc2d5dfSBarry Smith /*@
14058cc2d5dfSBarry Smith   PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
1406f1580f4eSBarry Smith   of multigrid when `PCMGType` is `PC_MG_MULTIPLICATIVE`
14078cc2d5dfSBarry Smith 
1408c3339decSBarry Smith   Logically Collective
14098cc2d5dfSBarry Smith 
14108cc2d5dfSBarry Smith   Input Parameters:
14118cc2d5dfSBarry Smith + pc - the multigrid context
14128cc2d5dfSBarry Smith - n  - number of cycles (default is 1)
14138cc2d5dfSBarry Smith 
14148cc2d5dfSBarry Smith   Options Database Key:
1415147403d9SBarry Smith . -pc_mg_multiplicative_cycles n - set the number of cycles
14168cc2d5dfSBarry Smith 
14178cc2d5dfSBarry Smith   Level: advanced
14188cc2d5dfSBarry Smith 
1419f1580f4eSBarry Smith   Note:
1420f1580f4eSBarry Smith   This is not associated with setting a v or w cycle, that is set with `PCMGSetCycleType()`
14218cc2d5dfSBarry Smith 
1422562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGSetCycleTypeOnLevel()`, `PCMGSetCycleType()`, `PCMGCycleType`, `PCMGType`
14238cc2d5dfSBarry Smith @*/
1424d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGMultiplicativeSetCycles(PC pc, PetscInt n)
1425d71ae5a4SJacob Faibussowitsch {
1426f3fbd535SBarry Smith   PC_MG *mg = (PC_MG *)pc->data;
14278cc2d5dfSBarry Smith 
14288cc2d5dfSBarry Smith   PetscFunctionBegin;
14290700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1430c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
14312a21e185SBarry Smith   mg->cyclesperpcapply = n;
14323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14338cc2d5dfSBarry Smith }
14348cc2d5dfSBarry Smith 
143566976f2fSJacob Faibussowitsch static PetscErrorCode PCMGSetGalerkin_MG(PC pc, PCMGGalerkinType use)
1436d71ae5a4SJacob Faibussowitsch {
1437fb15c04eSBarry Smith   PC_MG *mg = (PC_MG *)pc->data;
1438fb15c04eSBarry Smith 
1439fb15c04eSBarry Smith   PetscFunctionBegin;
14402134b1e4SBarry Smith   mg->galerkin = use;
14413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1442fb15c04eSBarry Smith }
1443fb15c04eSBarry Smith 
1444c2be2410SBarry Smith /*@
144597177400SBarry Smith   PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
144682b87d87SMatthew G. Knepley   finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1447c2be2410SBarry Smith 
1448c3339decSBarry Smith   Logically Collective
1449c2be2410SBarry Smith 
1450c2be2410SBarry Smith   Input Parameters:
1451c91913d3SJed Brown + pc  - the multigrid context
1452f1580f4eSBarry Smith - use - one of `PC_MG_GALERKIN_BOTH`, `PC_MG_GALERKIN_PMAT`, `PC_MG_GALERKIN_MAT`, or `PC_MG_GALERKIN_NONE`
1453c2be2410SBarry Smith 
1454c2be2410SBarry Smith   Options Database Key:
1455147403d9SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - set the matrices to form via the Galerkin process
1456c2be2410SBarry Smith 
1457c2be2410SBarry Smith   Level: intermediate
1458c2be2410SBarry Smith 
1459f1580f4eSBarry Smith   Note:
1460f1580f4eSBarry Smith   Some codes that use `PCMG` such as `PCGAMG` use Galerkin internally while constructing the hierarchy and thus do not
1461f1580f4eSBarry Smith   use the `PCMG` construction of the coarser grids.
14621c1aac46SBarry Smith 
1463562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetGalerkin()`, `PCMGGalerkinType`
1464c2be2410SBarry Smith @*/
1465d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetGalerkin(PC pc, PCMGGalerkinType use)
1466d71ae5a4SJacob Faibussowitsch {
1467c2be2410SBarry Smith   PetscFunctionBegin;
14680700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1469cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGSetGalerkin_C", (PC, PCMGGalerkinType), (pc, use));
14703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1471c2be2410SBarry Smith }
1472c2be2410SBarry Smith 
14733fc8bf9cSBarry Smith /*@
1474f1580f4eSBarry Smith   PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. A_i-1 = r_i * A_i * p_i
14753fc8bf9cSBarry Smith 
14763fc8bf9cSBarry Smith   Not Collective
14773fc8bf9cSBarry Smith 
14783fc8bf9cSBarry Smith   Input Parameter:
14793fc8bf9cSBarry Smith . pc - the multigrid context
14803fc8bf9cSBarry Smith 
14813fc8bf9cSBarry Smith   Output Parameter:
1482f1580f4eSBarry 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`
14833fc8bf9cSBarry Smith 
14843fc8bf9cSBarry Smith   Level: intermediate
14853fc8bf9cSBarry Smith 
1486562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetGalerkin()`, `PCMGGalerkinType`
14873fc8bf9cSBarry Smith @*/
1488d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetGalerkin(PC pc, PCMGGalerkinType *galerkin)
1489d71ae5a4SJacob Faibussowitsch {
1490f3fbd535SBarry Smith   PC_MG *mg = (PC_MG *)pc->data;
14913fc8bf9cSBarry Smith 
14923fc8bf9cSBarry Smith   PetscFunctionBegin;
14930700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
14942134b1e4SBarry Smith   *galerkin = mg->galerkin;
14953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14963fc8bf9cSBarry Smith }
14973fc8bf9cSBarry Smith 
149866976f2fSJacob Faibussowitsch static PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1499d71ae5a4SJacob Faibussowitsch {
1500f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *)pc->data;
1501f3b08a26SMatthew G. Knepley 
1502f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1503f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = adapt;
15043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1505f3b08a26SMatthew G. Knepley }
1506f3b08a26SMatthew G. Knepley 
150766976f2fSJacob Faibussowitsch static PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1508d71ae5a4SJacob Faibussowitsch {
1509f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *)pc->data;
1510f3b08a26SMatthew G. Knepley 
1511f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1512f3b08a26SMatthew G. Knepley   *adapt = mg->adaptInterpolation;
15133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1514f3b08a26SMatthew G. Knepley }
1515f3b08a26SMatthew G. Knepley 
151666976f2fSJacob Faibussowitsch static PetscErrorCode PCMGSetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType ctype)
1517d71ae5a4SJacob Faibussowitsch {
15182b3cbbdaSStefano Zampini   PC_MG *mg = (PC_MG *)pc->data;
15192b3cbbdaSStefano Zampini 
15202b3cbbdaSStefano Zampini   PetscFunctionBegin;
15212b3cbbdaSStefano Zampini   mg->adaptInterpolation = ctype != PCMG_ADAPT_NONE ? PETSC_TRUE : PETSC_FALSE;
15222b3cbbdaSStefano Zampini   mg->coarseSpaceType    = ctype;
15233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15242b3cbbdaSStefano Zampini }
15252b3cbbdaSStefano Zampini 
152666976f2fSJacob Faibussowitsch static PetscErrorCode PCMGGetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType *ctype)
1527d71ae5a4SJacob Faibussowitsch {
15282b3cbbdaSStefano Zampini   PC_MG *mg = (PC_MG *)pc->data;
15292b3cbbdaSStefano Zampini 
15302b3cbbdaSStefano Zampini   PetscFunctionBegin;
15312b3cbbdaSStefano Zampini   *ctype = mg->coarseSpaceType;
15323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15332b3cbbdaSStefano Zampini }
15342b3cbbdaSStefano Zampini 
153566976f2fSJacob Faibussowitsch static PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr)
1536d71ae5a4SJacob Faibussowitsch {
153741b6fd38SMatthew G. Knepley   PC_MG *mg = (PC_MG *)pc->data;
153841b6fd38SMatthew G. Knepley 
153941b6fd38SMatthew G. Knepley   PetscFunctionBegin;
154041b6fd38SMatthew G. Knepley   mg->compatibleRelaxation = cr;
15413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
154241b6fd38SMatthew G. Knepley }
154341b6fd38SMatthew G. Knepley 
154466976f2fSJacob Faibussowitsch static PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr)
1545d71ae5a4SJacob Faibussowitsch {
154641b6fd38SMatthew G. Knepley   PC_MG *mg = (PC_MG *)pc->data;
154741b6fd38SMatthew G. Knepley 
154841b6fd38SMatthew G. Knepley   PetscFunctionBegin;
154941b6fd38SMatthew G. Knepley   *cr = mg->compatibleRelaxation;
15503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
155141b6fd38SMatthew G. Knepley }
155241b6fd38SMatthew G. Knepley 
15532b3cbbdaSStefano Zampini /*@C
15542b3cbbdaSStefano Zampini   PCMGSetAdaptCoarseSpaceType - Set the type of adaptive coarse space.
15552b3cbbdaSStefano Zampini 
15562b3cbbdaSStefano 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.
15572b3cbbdaSStefano Zampini 
1558c3339decSBarry Smith   Logically Collective
15592b3cbbdaSStefano Zampini 
15602b3cbbdaSStefano Zampini   Input Parameters:
15612b3cbbdaSStefano Zampini + pc    - the multigrid context
15622b3cbbdaSStefano Zampini - ctype - the type of coarse space
15632b3cbbdaSStefano Zampini 
15642b3cbbdaSStefano Zampini   Options Database Keys:
15652b3cbbdaSStefano Zampini + -pc_mg_adapt_interp_n <int>             - The number of modes to use
15662b3cbbdaSStefano Zampini - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: none, polynomial, harmonic, eigenvector, generalized_eigenvector, gdsw
15672b3cbbdaSStefano Zampini 
15682b3cbbdaSStefano Zampini   Level: intermediate
15692b3cbbdaSStefano Zampini 
1570562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGCoarseSpaceType`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()`
15712b3cbbdaSStefano Zampini @*/
1572d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType ctype)
1573d71ae5a4SJacob Faibussowitsch {
15742b3cbbdaSStefano Zampini   PetscFunctionBegin;
15752b3cbbdaSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15762b3cbbdaSStefano Zampini   PetscValidLogicalCollectiveEnum(pc, ctype, 2);
15772b3cbbdaSStefano Zampini   PetscTryMethod(pc, "PCMGSetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType), (pc, ctype));
15783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15792b3cbbdaSStefano Zampini }
15802b3cbbdaSStefano Zampini 
15812b3cbbdaSStefano Zampini /*@C
15822b3cbbdaSStefano Zampini   PCMGGetAdaptCoarseSpaceType - Get the type of adaptive coarse space.
15832b3cbbdaSStefano Zampini 
15842b3cbbdaSStefano Zampini   Not Collective
15852b3cbbdaSStefano Zampini 
15862b3cbbdaSStefano Zampini   Input Parameter:
15872b3cbbdaSStefano Zampini . pc - the multigrid context
15882b3cbbdaSStefano Zampini 
15892b3cbbdaSStefano Zampini   Output Parameter:
15902b3cbbdaSStefano Zampini . ctype - the type of coarse space
15912b3cbbdaSStefano Zampini 
15922b3cbbdaSStefano Zampini   Level: intermediate
15932b3cbbdaSStefano Zampini 
1594562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGCoarseSpaceType`, `PCMGSetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()`
15952b3cbbdaSStefano Zampini @*/
1596d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType *ctype)
1597d71ae5a4SJacob Faibussowitsch {
15982b3cbbdaSStefano Zampini   PetscFunctionBegin;
15992b3cbbdaSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16004f572ea9SToby Isaac   PetscAssertPointer(ctype, 2);
16012b3cbbdaSStefano Zampini   PetscUseMethod(pc, "PCMGGetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType *), (pc, ctype));
16023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16032b3cbbdaSStefano Zampini }
16042b3cbbdaSStefano Zampini 
16052b3cbbdaSStefano Zampini /* MATT: REMOVE? */
1606f3b08a26SMatthew G. Knepley /*@
1607f3b08a26SMatthew 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.
1608f3b08a26SMatthew G. Knepley 
1609c3339decSBarry Smith   Logically Collective
1610f3b08a26SMatthew G. Knepley 
1611f3b08a26SMatthew G. Knepley   Input Parameters:
1612f3b08a26SMatthew G. Knepley + pc    - the multigrid context
1613f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator
1614f3b08a26SMatthew G. Knepley 
1615f3b08a26SMatthew G. Knepley   Options Database Keys:
1616f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp                     - Turn on adaptation
1617f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int>             - The number of modes to use, should be divisible by dimension
1618f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector
1619f3b08a26SMatthew G. Knepley 
1620f3b08a26SMatthew G. Knepley   Level: intermediate
1621f3b08a26SMatthew G. Knepley 
1622562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
1623f3b08a26SMatthew G. Knepley @*/
1624d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1625d71ae5a4SJacob Faibussowitsch {
1626f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1627f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1628cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGSetAdaptInterpolation_C", (PC, PetscBool), (pc, adapt));
16293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1630f3b08a26SMatthew G. Knepley }
1631f3b08a26SMatthew G. Knepley 
1632f3b08a26SMatthew G. Knepley /*@
1633f1580f4eSBarry Smith   PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh,
1634f1580f4eSBarry Smith   and thus accurately interpolated.
1635f3b08a26SMatthew G. Knepley 
163620f4b53cSBarry Smith   Not Collective
1637f3b08a26SMatthew G. Knepley 
1638f3b08a26SMatthew G. Knepley   Input Parameter:
1639f3b08a26SMatthew G. Knepley . pc - the multigrid context
1640f3b08a26SMatthew G. Knepley 
1641f3b08a26SMatthew G. Knepley   Output Parameter:
1642f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator
1643f3b08a26SMatthew G. Knepley 
1644f3b08a26SMatthew G. Knepley   Level: intermediate
1645f3b08a26SMatthew G. Knepley 
1646562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
1647f3b08a26SMatthew G. Knepley @*/
1648d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1649d71ae5a4SJacob Faibussowitsch {
1650f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1651f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16524f572ea9SToby Isaac   PetscAssertPointer(adapt, 2);
1653cac4c232SBarry Smith   PetscUseMethod(pc, "PCMGGetAdaptInterpolation_C", (PC, PetscBool *), (pc, adapt));
16543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1655f3b08a26SMatthew G. Knepley }
1656f3b08a26SMatthew G. Knepley 
16574b9ad928SBarry Smith /*@
165841b6fd38SMatthew G. Knepley   PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation.
165941b6fd38SMatthew G. Knepley 
1660c3339decSBarry Smith   Logically Collective
166141b6fd38SMatthew G. Knepley 
166241b6fd38SMatthew G. Knepley   Input Parameters:
166341b6fd38SMatthew G. Knepley + pc - the multigrid context
166441b6fd38SMatthew G. Knepley - cr - flag for compatible relaxation
166541b6fd38SMatthew G. Knepley 
1666f1580f4eSBarry Smith   Options Database Key:
166741b6fd38SMatthew G. Knepley . -pc_mg_adapt_cr - Turn on compatible relaxation
166841b6fd38SMatthew G. Knepley 
166941b6fd38SMatthew G. Knepley   Level: intermediate
167041b6fd38SMatthew G. Knepley 
1671562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetAdaptCR()`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
167241b6fd38SMatthew G. Knepley @*/
1673d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr)
1674d71ae5a4SJacob Faibussowitsch {
167541b6fd38SMatthew G. Knepley   PetscFunctionBegin;
167641b6fd38SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1677cac4c232SBarry Smith   PetscTryMethod(pc, "PCMGSetAdaptCR_C", (PC, PetscBool), (pc, cr));
16783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
167941b6fd38SMatthew G. Knepley }
168041b6fd38SMatthew G. Knepley 
168141b6fd38SMatthew G. Knepley /*@
168241b6fd38SMatthew G. Knepley   PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation.
168341b6fd38SMatthew G. Knepley 
168420f4b53cSBarry Smith   Not Collective
168541b6fd38SMatthew G. Knepley 
168641b6fd38SMatthew G. Knepley   Input Parameter:
168741b6fd38SMatthew G. Knepley . pc - the multigrid context
168841b6fd38SMatthew G. Knepley 
168941b6fd38SMatthew G. Knepley   Output Parameter:
169041b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion
169141b6fd38SMatthew G. Knepley 
169241b6fd38SMatthew G. Knepley   Level: intermediate
169341b6fd38SMatthew G. Knepley 
1694562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
169541b6fd38SMatthew G. Knepley @*/
1696d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr)
1697d71ae5a4SJacob Faibussowitsch {
169841b6fd38SMatthew G. Knepley   PetscFunctionBegin;
169941b6fd38SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
17004f572ea9SToby Isaac   PetscAssertPointer(cr, 2);
1701cac4c232SBarry Smith   PetscUseMethod(pc, "PCMGGetAdaptCR_C", (PC, PetscBool *), (pc, cr));
17023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
170341b6fd38SMatthew G. Knepley }
170441b6fd38SMatthew G. Knepley 
170541b6fd38SMatthew G. Knepley /*@
170606792cafSBarry Smith   PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1707f1580f4eSBarry Smith   on all levels.  Use `PCMGDistinctSmoothUp()` to create separate up and down smoothers if you want different numbers of
1708710315b6SLawrence Mitchell   pre- and post-smoothing steps.
170906792cafSBarry Smith 
1710c3339decSBarry Smith   Logically Collective
171106792cafSBarry Smith 
171206792cafSBarry Smith   Input Parameters:
1713feefa0e1SJacob Faibussowitsch + pc - the multigrid context
171406792cafSBarry Smith - n  - the number of smoothing steps
171506792cafSBarry Smith 
171606792cafSBarry Smith   Options Database Key:
1717a2b725a8SWilliam Gropp . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
171806792cafSBarry Smith 
171906792cafSBarry Smith   Level: advanced
172006792cafSBarry Smith 
1721f1580f4eSBarry Smith   Note:
1722f1580f4eSBarry 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.
172306792cafSBarry Smith 
1724562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetDistinctSmoothUp()`
172506792cafSBarry Smith @*/
1726d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetNumberSmooth(PC pc, PetscInt n)
1727d71ae5a4SJacob Faibussowitsch {
172806792cafSBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
172906792cafSBarry Smith   PC_MG_Levels **mglevels = mg->levels;
173006792cafSBarry Smith   PetscInt       i, levels;
173106792cafSBarry Smith 
173206792cafSBarry Smith   PetscFunctionBegin;
173306792cafSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
173406792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc, n, 2);
173528b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling");
173606792cafSBarry Smith   levels = mglevels[0]->levels;
173706792cafSBarry Smith 
173806792cafSBarry Smith   for (i = 1; i < levels; i++) {
17399566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[i]->smoothu, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, n));
17409566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, n));
174106792cafSBarry Smith     mg->default_smoothu = n;
174206792cafSBarry Smith     mg->default_smoothd = n;
174306792cafSBarry Smith   }
17443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
174506792cafSBarry Smith }
174606792cafSBarry Smith 
1747f442ab6aSBarry Smith /*@
1748f1580f4eSBarry Smith   PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate `KSP` from the down (pre) smoother on all levels
1749710315b6SLawrence Mitchell   and adds the suffix _up to the options name
1750f442ab6aSBarry Smith 
1751c3339decSBarry Smith   Logically Collective
1752f442ab6aSBarry Smith 
1753f1580f4eSBarry Smith   Input Parameter:
1754f442ab6aSBarry Smith . pc - the preconditioner context
1755f442ab6aSBarry Smith 
1756f442ab6aSBarry Smith   Options Database Key:
1757147403d9SBarry Smith . -pc_mg_distinct_smoothup <bool> - use distinct smoothing objects
1758f442ab6aSBarry Smith 
1759f442ab6aSBarry Smith   Level: advanced
1760f442ab6aSBarry Smith 
1761f1580f4eSBarry Smith   Note:
1762f1580f4eSBarry 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.
1763f442ab6aSBarry Smith 
1764562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGSetNumberSmooth()`
1765f442ab6aSBarry Smith @*/
1766d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGSetDistinctSmoothUp(PC pc)
1767d71ae5a4SJacob Faibussowitsch {
1768f442ab6aSBarry Smith   PC_MG         *mg       = (PC_MG *)pc->data;
1769f442ab6aSBarry Smith   PC_MG_Levels **mglevels = mg->levels;
1770f442ab6aSBarry Smith   PetscInt       i, levels;
1771f442ab6aSBarry Smith   KSP            subksp;
1772f442ab6aSBarry Smith 
1773f442ab6aSBarry Smith   PetscFunctionBegin;
1774f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
177528b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling");
1776f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1777f442ab6aSBarry Smith 
1778f442ab6aSBarry Smith   for (i = 1; i < levels; i++) {
1779710315b6SLawrence Mitchell     const char *prefix = NULL;
1780f442ab6aSBarry Smith     /* make sure smoother up and down are different */
17819566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmootherUp(pc, i, &subksp));
17829566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(mglevels[i]->smoothd, &prefix));
17839566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(subksp, prefix));
17849566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(subksp, "up_"));
1785f442ab6aSBarry Smith   }
17863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1787f442ab6aSBarry Smith }
1788f442ab6aSBarry Smith 
178907a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
179066976f2fSJacob Faibussowitsch static PetscErrorCode PCGetInterpolations_MG(PC pc, PetscInt *num_levels, Mat *interpolations[])
1791d71ae5a4SJacob Faibussowitsch {
179207a4832bSFande Kong   PC_MG         *mg       = (PC_MG *)pc->data;
179307a4832bSFande Kong   PC_MG_Levels **mglevels = mg->levels;
179407a4832bSFande Kong   Mat           *mat;
179507a4832bSFande Kong   PetscInt       l;
179607a4832bSFande Kong 
179707a4832bSFande Kong   PetscFunctionBegin;
179828b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
17999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mg->nlevels, &mat));
180007a4832bSFande Kong   for (l = 1; l < mg->nlevels; l++) {
180107a4832bSFande Kong     mat[l - 1] = mglevels[l]->interpolate;
18029566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat[l - 1]));
180307a4832bSFande Kong   }
180407a4832bSFande Kong   *num_levels     = mg->nlevels;
180507a4832bSFande Kong   *interpolations = mat;
18063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
180707a4832bSFande Kong }
180807a4832bSFande Kong 
180907a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
181066976f2fSJacob Faibussowitsch static PetscErrorCode PCGetCoarseOperators_MG(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
1811d71ae5a4SJacob Faibussowitsch {
181207a4832bSFande Kong   PC_MG         *mg       = (PC_MG *)pc->data;
181307a4832bSFande Kong   PC_MG_Levels **mglevels = mg->levels;
181407a4832bSFande Kong   PetscInt       l;
181507a4832bSFande Kong   Mat           *mat;
181607a4832bSFande Kong 
181707a4832bSFande Kong   PetscFunctionBegin;
181828b400f6SJacob Faibussowitsch   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
18199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mg->nlevels, &mat));
182007a4832bSFande Kong   for (l = 0; l < mg->nlevels - 1; l++) {
1821*f4f49eeaSPierre Jolivet     PetscCall(KSPGetOperators(mglevels[l]->smoothd, NULL, &mat[l]));
18229566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat[l]));
182307a4832bSFande Kong   }
182407a4832bSFande Kong   *num_levels      = mg->nlevels;
182507a4832bSFande Kong   *coarseOperators = mat;
18263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
182707a4832bSFande Kong }
182807a4832bSFande Kong 
1829f3b08a26SMatthew G. Knepley /*@C
1830f1580f4eSBarry Smith   PCMGRegisterCoarseSpaceConstructor -  Adds a method to the `PCMG` package for coarse space construction.
1831f3b08a26SMatthew G. Knepley 
183220f4b53cSBarry Smith   Not Collective
1833f3b08a26SMatthew G. Knepley 
1834f3b08a26SMatthew G. Knepley   Input Parameters:
1835f3b08a26SMatthew G. Knepley + name     - name of the constructor
1836f3b08a26SMatthew G. Knepley - function - constructor routine
1837f3b08a26SMatthew G. Knepley 
183820f4b53cSBarry Smith   Calling sequence of `function`:
183920f4b53cSBarry Smith + pc        - The `PC` object
1840f1580f4eSBarry Smith . l         - The multigrid level, 0 is the coarse level
184120f4b53cSBarry Smith . dm        - The `DM` for this level
1842f1580f4eSBarry Smith . smooth    - The level smoother
1843f1580f4eSBarry Smith . Nc        - The size of the coarse space
1844f1580f4eSBarry Smith . initGuess - Basis for an initial guess for the space
1845f1580f4eSBarry Smith - coarseSp  - A basis for the computed coarse space
1846f3b08a26SMatthew G. Knepley 
1847f3b08a26SMatthew G. Knepley   Level: advanced
1848f3b08a26SMatthew G. Knepley 
1849feefa0e1SJacob Faibussowitsch   Developer Notes:
1850f1580f4eSBarry Smith   How come this is not used by `PCGAMG`?
1851f1580f4eSBarry Smith 
1852562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGGetCoarseSpaceConstructor()`, `PCRegister()`
1853f3b08a26SMatthew G. Knepley @*/
185404c3f3b8SBarry Smith PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp))
1855d71ae5a4SJacob Faibussowitsch {
1856f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
18579566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
18589566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCMGCoarseList, name, function));
18593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1860f3b08a26SMatthew G. Knepley }
1861f3b08a26SMatthew G. Knepley 
1862f3b08a26SMatthew G. Knepley /*@C
1863f3b08a26SMatthew G. Knepley   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.
1864f3b08a26SMatthew G. Knepley 
186520f4b53cSBarry Smith   Not Collective
1866f3b08a26SMatthew G. Knepley 
1867f3b08a26SMatthew G. Knepley   Input Parameter:
1868f3b08a26SMatthew G. Knepley . name - name of the constructor
1869f3b08a26SMatthew G. Knepley 
1870f3b08a26SMatthew G. Knepley   Output Parameter:
1871f3b08a26SMatthew G. Knepley . function - constructor routine
1872f3b08a26SMatthew G. Knepley 
1873f3b08a26SMatthew G. Knepley   Level: advanced
1874f3b08a26SMatthew G. Knepley 
1875562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMG`, `PCMGRegisterCoarseSpaceConstructor()`, `PCRegister()`
1876f3b08a26SMatthew G. Knepley @*/
1877d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *))
1878d71ae5a4SJacob Faibussowitsch {
1879f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
18809566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PCMGCoarseList, name, function));
18813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1882f3b08a26SMatthew G. Knepley }
1883f3b08a26SMatthew G. Knepley 
18843b09bd56SBarry Smith /*MC
1885ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
18863b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
18873b09bd56SBarry Smith 
18883b09bd56SBarry Smith    Options Database Keys:
18893b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1890391689abSStefano Zampini .  -pc_mg_cycle_type <v,w> - provide the cycle desired
18918c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
18923b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1893710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
18942134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
18958cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
18968cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1897e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
18988cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
18998cf6037eSBarry Smith                         to the binary output file called binaryoutput
19003b09bd56SBarry Smith 
190120f4b53cSBarry Smith    Level: intermediate
190220f4b53cSBarry Smith 
190395452b02SPatrick Sanan    Notes:
190404c3f3b8SBarry Smith    The Krylov solver (if any) and preconditioner (smoother) and their parameters are controlled from the options database with the standard
19058f4fb22eSMark Adams    options database keywords prefixed with `-mg_levels_` to affect all the levels but the coarsest, which is controlled with `-mg_coarse_`,
19068f4fb22eSMark Adams    and the finest where `-mg_fine_` can override `-mg_levels_`.  One can set different preconditioners etc on specific levels with the prefix
19078f4fb22eSMark Adams    `-mg_levels_n_` where `n` is the level number (zero being the coarse level. For example
190804c3f3b8SBarry Smith .vb
190904c3f3b8SBarry Smith    -mg_levels_ksp_type gmres -mg_levels_pc_type bjacobi -mg_coarse_pc_type svd -mg_levels_2_pc_type sor
191004c3f3b8SBarry Smith .ve
191104c3f3b8SBarry Smith    These options also work for controlling the smoothers etc inside `PCGAMG`
191204c3f3b8SBarry Smith 
1913f1580f4eSBarry 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
19143b09bd56SBarry Smith 
19158cf6037eSBarry Smith    When run with a single level the smoother options are used on that level NOT the coarse grid solver options
19168cf6037eSBarry Smith 
1917f1580f4eSBarry Smith    When run with `KSPRICHARDSON` the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
191823067569SBarry Smith    is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
191923067569SBarry Smith    (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
192023067569SBarry Smith    residual is computed at the end of each cycle.
192123067569SBarry Smith 
192204c3f3b8SBarry Smith .seealso: [](sec_mg), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCEXOTIC`, `PCGAMG`, `PCML`, `PCHYPRE`
1923db781477SPatrick Sanan           `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `PCMGSetCycleType()`,
1924db781477SPatrick Sanan           `PCMGSetDistinctSmoothUp()`, `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`,
1925db781477SPatrick Sanan           `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`,
1926f1580f4eSBarry Smith           `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`,
1927f1580f4eSBarry Smith           `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
19283b09bd56SBarry Smith M*/
19293b09bd56SBarry Smith 
1930d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1931d71ae5a4SJacob Faibussowitsch {
1932f3fbd535SBarry Smith   PC_MG *mg;
1933f3fbd535SBarry Smith 
19344b9ad928SBarry Smith   PetscFunctionBegin;
19354dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&mg));
19363ec1f749SStefano Zampini   pc->data               = mg;
1937f3fbd535SBarry Smith   mg->nlevels            = -1;
193810eca3edSLisandro Dalcin   mg->am                 = PC_MG_MULTIPLICATIVE;
19392134b1e4SBarry Smith   mg->galerkin           = PC_MG_GALERKIN_NONE;
1940f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = PETSC_FALSE;
1941f3b08a26SMatthew G. Knepley   mg->Nc                 = -1;
1942f3b08a26SMatthew G. Knepley   mg->eigenvalue         = -1;
1943f3fbd535SBarry Smith 
194437a44384SMark Adams   pc->useAmat = PETSC_TRUE;
194537a44384SMark Adams 
19464b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
1947fcb023d4SJed Brown   pc->ops->applytranspose = PCApplyTranspose_MG;
194830b0564aSStefano Zampini   pc->ops->matapply       = PCMatApply_MG;
19494b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1950a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
19514b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
19524b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
19534b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1954fb15c04eSBarry Smith 
19559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataRegister(&mg->eigenvalue));
19569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", PCMGSetGalerkin_MG));
19579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", PCMGGetLevels_MG));
19589566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", PCMGSetLevels_MG));
19599566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", PCGetInterpolations_MG));
19609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", PCGetCoarseOperators_MG));
19619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", PCMGSetAdaptInterpolation_MG));
19629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", PCMGGetAdaptInterpolation_MG));
19639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", PCMGSetAdaptCR_MG));
19649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", PCMGGetAdaptCR_MG));
19652b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", PCMGSetAdaptCoarseSpaceType_MG));
19662b3cbbdaSStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", PCMGGetAdaptCoarseSpaceType_MG));
19673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19684b9ad928SBarry Smith }
1969