xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 7827d75ba736e00c19d105c058b1c2ddcca945c7)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>                    /*I "petscksp.h" I*/
641b6fd38SMatthew G. Knepley #include <petsc/private/kspimpl.h>
71e25c274SJed Brown #include <petscdm.h>
8391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*);
94b9ad928SBarry Smith 
10f3b08a26SMatthew G. Knepley /*
11f3b08a26SMatthew G. Knepley    Contains the list of registered coarse space construction routines
12f3b08a26SMatthew G. Knepley */
13f3b08a26SMatthew G. Knepley PetscFunctionList PCMGCoarseList = NULL;
14f3b08a26SMatthew G. Knepley 
1530b0564aSStefano Zampini PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PetscBool transpose,PetscBool matapp,PCRichardsonConvergedReason *reason)
164b9ad928SBarry Smith {
1731567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1831567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
1931567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
204b9ad928SBarry Smith 
214b9ad928SBarry Smith   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0));
23fcb023d4SJed Brown   if (!transpose) {
2430b0564aSStefano Zampini     if (matapp) {
259566063dSJacob Faibussowitsch       PetscCall(KSPMatSolve(mglevels->smoothd,mglevels->B,mglevels->X));  /* pre-smooth */
269566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->smoothd,pc,NULL));
2730b0564aSStefano Zampini     } else {
289566063dSJacob Faibussowitsch       PetscCall(KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x));  /* pre-smooth */
299566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->smoothd,pc,mglevels->x));
3030b0564aSStefano Zampini     }
31fcb023d4SJed Brown   } else {
3228b400f6SJacob Faibussowitsch     PetscCheck(!matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported");
339566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(mglevels->smoothu,mglevels->b,mglevels->x)); /* transpose of post-smooth */
349566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(mglevels->smoothu,pc,mglevels->x));
35fcb023d4SJed Brown   }
369566063dSJacob Faibussowitsch   if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0));
3731567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
389566063dSJacob Faibussowitsch     if (mglevels->eventresidual) PetscCall(PetscLogEventBegin(mglevels->eventresidual,0,0,0,0));
3930b0564aSStefano Zampini     if (matapp && !mglevels->R) {
409566063dSJacob Faibussowitsch       PetscCall(MatDuplicate(mglevels->B,MAT_DO_NOT_COPY_VALUES,&mglevels->R));
4130b0564aSStefano Zampini     }
42fcb023d4SJed Brown     if (!transpose) {
439566063dSJacob Faibussowitsch       if (matapp) PetscCall((*mglevels->matresidual)(mglevels->A,mglevels->B,mglevels->X,mglevels->R));
449566063dSJacob Faibussowitsch       else PetscCall((*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r));
45fcb023d4SJed Brown     } else {
469566063dSJacob Faibussowitsch       if (matapp) PetscCall((*mglevels->matresidualtranspose)(mglevels->A,mglevels->B,mglevels->X,mglevels->R));
479566063dSJacob Faibussowitsch       else PetscCall((*mglevels->residualtranspose)(mglevels->A,mglevels->b,mglevels->x,mglevels->r));
48fcb023d4SJed Brown     }
499566063dSJacob Faibussowitsch     if (mglevels->eventresidual) PetscCall(PetscLogEventEnd(mglevels->eventresidual,0,0,0,0));
504b9ad928SBarry Smith 
514b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
5231567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
534b9ad928SBarry Smith       PetscReal rnorm;
549566063dSJacob Faibussowitsch       PetscCall(VecNorm(mglevels->r,NORM_2,&rnorm));
554b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
5670441072SBarry Smith         if (rnorm < mg->abstol) {
574d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
589566063dSJacob Faibussowitsch           PetscCall(PetscInfo(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol));
594b9ad928SBarry Smith         } else {
604d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
619566063dSJacob Faibussowitsch           PetscCall(PetscInfo(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol));
624b9ad928SBarry Smith         }
634b9ad928SBarry Smith         PetscFunctionReturn(0);
644b9ad928SBarry Smith       }
654b9ad928SBarry Smith     }
664b9ad928SBarry Smith 
6731567311SBarry Smith     mgc = *(mglevelsin - 1);
689566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0));
69fcb023d4SJed Brown     if (!transpose) {
709566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatRestrict(mglevels->restrct,mglevels->R,&mgc->B));
719566063dSJacob Faibussowitsch       else PetscCall(MatRestrict(mglevels->restrct,mglevels->r,mgc->b));
72fcb023d4SJed Brown     } else {
739566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatRestrict(mglevels->interpolate,mglevels->R,&mgc->B));
749566063dSJacob Faibussowitsch       else PetscCall(MatRestrict(mglevels->interpolate,mglevels->r,mgc->b));
75fcb023d4SJed Brown     }
769566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0));
7730b0564aSStefano Zampini     if (matapp) {
7830b0564aSStefano Zampini       if (!mgc->X) {
799566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(mgc->B,MAT_DO_NOT_COPY_VALUES,&mgc->X));
8030b0564aSStefano Zampini       } else {
819566063dSJacob Faibussowitsch         PetscCall(MatZeroEntries(mgc->X));
8230b0564aSStefano Zampini       }
8330b0564aSStefano Zampini     } else {
849566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(mgc->x));
8530b0564aSStefano Zampini     }
864b9ad928SBarry Smith     while (cycles--) {
879566063dSJacob Faibussowitsch       PetscCall(PCMGMCycle_Private(pc,mglevelsin-1,transpose,matapp,reason));
884b9ad928SBarry Smith     }
899566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0));
90fcb023d4SJed Brown     if (!transpose) {
919566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatInterpolateAdd(mglevels->interpolate,mgc->X,mglevels->X,&mglevels->X));
929566063dSJacob Faibussowitsch       else PetscCall(MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x));
93fcb023d4SJed Brown     } else {
949566063dSJacob Faibussowitsch       PetscCall(MatInterpolateAdd(mglevels->restrct,mgc->x,mglevels->x,mglevels->x));
95fcb023d4SJed Brown     }
969566063dSJacob Faibussowitsch     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0));
979566063dSJacob Faibussowitsch     if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0));
98fcb023d4SJed Brown     if (!transpose) {
9930b0564aSStefano Zampini       if (matapp) {
1009566063dSJacob Faibussowitsch         PetscCall(KSPMatSolve(mglevels->smoothu,mglevels->B,mglevels->X));    /* post smooth */
1019566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(mglevels->smoothu,pc,NULL));
10230b0564aSStefano Zampini       } else {
1039566063dSJacob Faibussowitsch         PetscCall(KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x));    /* post smooth */
1049566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(mglevels->smoothu,pc,mglevels->x));
10530b0564aSStefano Zampini       }
106fcb023d4SJed Brown     } else {
10728b400f6SJacob Faibussowitsch       PetscCheck(!matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported");
1089566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(mglevels->smoothd,mglevels->b,mglevels->x));    /* post smooth */
1099566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->smoothd,pc,mglevels->x));
110fcb023d4SJed Brown     }
11141b6fd38SMatthew G. Knepley     if (mglevels->cr) {
11228b400f6SJacob Faibussowitsch       PetscCheck(!matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported");
11341b6fd38SMatthew G. Knepley       /* TODO Turn on copy and turn off noisy if we have an exact solution
1149566063dSJacob Faibussowitsch       PetscCall(VecCopy(mglevels->x, mglevels->crx));
1159566063dSJacob Faibussowitsch       PetscCall(VecCopy(mglevels->b, mglevels->crb)); */
1169566063dSJacob Faibussowitsch       PetscCall(KSPSetNoisy_Private(mglevels->crx));
1179566063dSJacob Faibussowitsch       PetscCall(KSPSolve(mglevels->cr,mglevels->crb,mglevels->crx));    /* compatible relaxation */
1189566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels->cr,pc,mglevels->crx));
11941b6fd38SMatthew G. Knepley     }
1209566063dSJacob Faibussowitsch     if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0));
1214b9ad928SBarry Smith   }
1224b9ad928SBarry Smith   PetscFunctionReturn(0);
1234b9ad928SBarry Smith }
1244b9ad928SBarry Smith 
125ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
1264b9ad928SBarry Smith {
127f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
128f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
129391689abSStefano Zampini   PC             tpc;
130391689abSStefano Zampini   PetscBool      changeu,changed;
131f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
1324b9ad928SBarry Smith 
1334b9ad928SBarry Smith   PetscFunctionBegin;
134a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
135a762d673SBarry Smith   for (i=0; i<levels; i++) {
136a762d673SBarry Smith     if (!mglevels[i]->A) {
1379566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL));
1389566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A));
139a762d673SBarry Smith     }
140a762d673SBarry Smith   }
141391689abSStefano Zampini 
1429566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[levels-1]->smoothd,&tpc));
1439566063dSJacob Faibussowitsch   PetscCall(PCPreSolveChangeRHS(tpc,&changed));
1449566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[levels-1]->smoothu,&tpc));
1459566063dSJacob Faibussowitsch   PetscCall(PCPreSolveChangeRHS(tpc,&changeu));
146391689abSStefano Zampini   if (!changed && !changeu) {
1479566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[levels-1]->b));
148f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
149391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
150391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
151391689abSStefano Zampini       Vec *vec;
152391689abSStefano Zampini 
1539566063dSJacob Faibussowitsch       PetscCall(KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL));
154391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
1559566063dSJacob Faibussowitsch       PetscCall(PetscFree(vec));
156391689abSStefano Zampini     }
1579566063dSJacob Faibussowitsch     PetscCall(VecCopy(b,mglevels[levels-1]->b));
158391689abSStefano Zampini   }
159f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
1604b9ad928SBarry Smith 
16131567311SBarry Smith   mg->rtol   = rtol;
16231567311SBarry Smith   mg->abstol = abstol;
16331567311SBarry Smith   mg->dtol   = dtol;
1644b9ad928SBarry Smith   if (rtol) {
1654b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
1664b9ad928SBarry Smith     PetscReal rnorm;
1677319c654SBarry Smith     if (zeroguess) {
1689566063dSJacob Faibussowitsch       PetscCall(VecNorm(b,NORM_2,&rnorm));
1697319c654SBarry Smith     } else {
1709566063dSJacob Faibussowitsch       PetscCall((*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w));
1719566063dSJacob Faibussowitsch       PetscCall(VecNorm(w,NORM_2,&rnorm));
1727319c654SBarry Smith     }
17331567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
1742fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
1752fa5cd67SKarl Rupp   else mg->ttol = 0.0;
1764b9ad928SBarry Smith 
1774d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
178336babb1SJed Brown      stop prematurely due to small residual */
1794d0a8057SBarry Smith   for (i=1; i<levels; i++) {
1809566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
181f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
18223067569SBarry Smith       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
1839566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE));
1849566063dSJacob Faibussowitsch       PetscCall(KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
1854b9ad928SBarry Smith     }
1864d0a8057SBarry Smith   }
1874d0a8057SBarry Smith 
1884d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1894d0a8057SBarry Smith   for (i=0; i<its; i++) {
1909566063dSJacob Faibussowitsch     PetscCall(PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_FALSE,PETSC_FALSE,reason));
1914d0a8057SBarry Smith     if (*reason) break;
1924d0a8057SBarry Smith   }
1934d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1944d0a8057SBarry Smith   *outits = i;
195391689abSStefano Zampini   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
1964b9ad928SBarry Smith   PetscFunctionReturn(0);
1974b9ad928SBarry Smith }
1984b9ad928SBarry Smith 
1993aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
2003aeaf226SBarry Smith {
2013aeaf226SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
2023aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
203f3b08a26SMatthew G. Knepley   PetscInt       i,c,n;
2043aeaf226SBarry Smith 
2053aeaf226SBarry Smith   PetscFunctionBegin;
2063aeaf226SBarry Smith   if (mglevels) {
2073aeaf226SBarry Smith     n = mglevels[0]->levels;
2083aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
2099566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i+1]->r));
2109566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->b));
2119566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->x));
2129566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i+1]->R));
2139566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i]->B));
2149566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i]->X));
2159566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->crx));
2169566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i]->crb));
2179566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i+1]->restrct));
2189566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i+1]->interpolate));
2199566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i+1]->inject));
2209566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[i+1]->rscale));
2213aeaf226SBarry Smith     }
2229566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[n-1]->crx));
2239566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[n-1]->crb));
224391689abSStefano Zampini     /* this is not null only if the smoother on the finest level
225391689abSStefano Zampini        changes the rhs during PreSolve */
2269566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mglevels[n-1]->b));
2279566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&mglevels[n-1]->B));
2283aeaf226SBarry Smith 
2293aeaf226SBarry Smith     for (i=0; i<n; i++) {
2309566063dSJacob Faibussowitsch       if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) PetscCall(VecDestroy(&mglevels[i]->coarseSpace[c]));
2319566063dSJacob Faibussowitsch       PetscCall(PetscFree(mglevels[i]->coarseSpace));
232f3b08a26SMatthew G. Knepley       mglevels[i]->coarseSpace = NULL;
2339566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[i]->A));
2343aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2359566063dSJacob Faibussowitsch         PetscCall(KSPReset(mglevels[i]->smoothd));
2363aeaf226SBarry Smith       }
2379566063dSJacob Faibussowitsch       PetscCall(KSPReset(mglevels[i]->smoothu));
2389566063dSJacob Faibussowitsch       if (mglevels[i]->cr) PetscCall(KSPReset(mglevels[i]->cr));
2393aeaf226SBarry Smith     }
240f3b08a26SMatthew G. Knepley     mg->Nc = 0;
2413aeaf226SBarry Smith   }
2423aeaf226SBarry Smith   PetscFunctionReturn(0);
2433aeaf226SBarry Smith }
2443aeaf226SBarry Smith 
24541b6fd38SMatthew G. Knepley /* Implementing CR
24641b6fd38SMatthew G. Knepley 
24741b6fd38SMatthew 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
24841b6fd38SMatthew G. Knepley 
24941b6fd38SMatthew G. Knepley   Inj^T (Inj Inj^T)^{-1} Inj
25041b6fd38SMatthew G. Knepley 
25141b6fd38SMatthew G. Knepley and if Inj is a VecScatter, as it is now in PETSc, we have
25241b6fd38SMatthew G. Knepley 
25341b6fd38SMatthew G. Knepley   Inj^T Inj
25441b6fd38SMatthew G. Knepley 
25541b6fd38SMatthew G. Knepley and
25641b6fd38SMatthew G. Knepley 
25741b6fd38SMatthew G. Knepley   S = I - Inj^T Inj
25841b6fd38SMatthew G. Knepley 
25941b6fd38SMatthew G. Knepley since
26041b6fd38SMatthew G. Knepley 
26141b6fd38SMatthew G. Knepley   Inj S = Inj - (Inj Inj^T) Inj = 0.
26241b6fd38SMatthew G. Knepley 
26341b6fd38SMatthew G. Knepley Brannick suggests
26441b6fd38SMatthew G. Knepley 
26541b6fd38SMatthew G. Knepley   A \to S^T A S  \qquad\mathrm{and}\qquad M \to S^T M S
26641b6fd38SMatthew G. Knepley 
26741b6fd38SMatthew 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
26841b6fd38SMatthew G. Knepley 
26941b6fd38SMatthew G. Knepley   M^{-1} A \to S M^{-1} A S
27041b6fd38SMatthew G. Knepley 
27141b6fd38SMatthew G. Knepley In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left.
27241b6fd38SMatthew G. Knepley 
27341b6fd38SMatthew G. Knepley   Check: || Inj P - I ||_F < tol
27441b6fd38SMatthew G. Knepley   Check: In general, Inj Inj^T = I
27541b6fd38SMatthew G. Knepley */
27641b6fd38SMatthew G. Knepley 
27741b6fd38SMatthew G. Knepley typedef struct {
27841b6fd38SMatthew G. Knepley   PC       mg;  /* The PCMG object */
27941b6fd38SMatthew G. Knepley   PetscInt l;   /* The multigrid level for this solver */
28041b6fd38SMatthew G. Knepley   Mat      Inj; /* The injection matrix */
28141b6fd38SMatthew G. Knepley   Mat      S;   /* I - Inj^T Inj */
28241b6fd38SMatthew G. Knepley } CRContext;
28341b6fd38SMatthew G. Knepley 
28441b6fd38SMatthew G. Knepley static PetscErrorCode CRSetup_Private(PC pc)
28541b6fd38SMatthew G. Knepley {
28641b6fd38SMatthew G. Knepley   CRContext     *ctx;
28741b6fd38SMatthew G. Knepley   Mat            It;
28841b6fd38SMatthew G. Knepley 
28941b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
2909566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &ctx));
2919566063dSJacob Faibussowitsch   PetscCall(PCMGGetInjection(ctx->mg, ctx->l, &It));
29228b400f6SJacob Faibussowitsch   PetscCheck(It,PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG");
2939566063dSJacob Faibussowitsch   PetscCall(MatCreateTranspose(It, &ctx->Inj));
2949566063dSJacob Faibussowitsch   PetscCall(MatCreateNormal(ctx->Inj, &ctx->S));
2959566063dSJacob Faibussowitsch   PetscCall(MatScale(ctx->S, -1.0));
2969566063dSJacob Faibussowitsch   PetscCall(MatShift(ctx->S,  1.0));
29741b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
29841b6fd38SMatthew G. Knepley }
29941b6fd38SMatthew G. Knepley 
30041b6fd38SMatthew G. Knepley static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y)
30141b6fd38SMatthew G. Knepley {
30241b6fd38SMatthew G. Knepley   CRContext     *ctx;
30341b6fd38SMatthew G. Knepley 
30441b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
3059566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &ctx));
3069566063dSJacob Faibussowitsch   PetscCall(MatMult(ctx->S, x, y));
30741b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
30841b6fd38SMatthew G. Knepley }
30941b6fd38SMatthew G. Knepley 
31041b6fd38SMatthew G. Knepley static PetscErrorCode CRDestroy_Private(PC pc)
31141b6fd38SMatthew G. Knepley {
31241b6fd38SMatthew G. Knepley   CRContext     *ctx;
31341b6fd38SMatthew G. Knepley 
31441b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
3159566063dSJacob Faibussowitsch   PetscCall(PCShellGetContext(pc, &ctx));
3169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->Inj));
3179566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx->S));
3189566063dSJacob Faibussowitsch   PetscCall(PetscFree(ctx));
3199566063dSJacob Faibussowitsch   PetscCall(PCShellSetContext(pc, NULL));
32041b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
32141b6fd38SMatthew G. Knepley }
32241b6fd38SMatthew G. Knepley 
32341b6fd38SMatthew G. Knepley static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr)
32441b6fd38SMatthew G. Knepley {
32541b6fd38SMatthew G. Knepley   CRContext     *ctx;
32641b6fd38SMatthew G. Knepley 
32741b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
3289566063dSJacob Faibussowitsch   PetscCall(PCCreate(PetscObjectComm((PetscObject) pc), cr));
3299566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *cr, "S (complementary projector to injection)"));
3309566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(1, &ctx));
33141b6fd38SMatthew G. Knepley   ctx->mg = pc;
33241b6fd38SMatthew G. Knepley   ctx->l  = l;
3339566063dSJacob Faibussowitsch   PetscCall(PCSetType(*cr, PCSHELL));
3349566063dSJacob Faibussowitsch   PetscCall(PCShellSetContext(*cr, ctx));
3359566063dSJacob Faibussowitsch   PetscCall(PCShellSetApply(*cr, CRApply_Private));
3369566063dSJacob Faibussowitsch   PetscCall(PCShellSetSetUp(*cr, CRSetup_Private));
3379566063dSJacob Faibussowitsch   PetscCall(PCShellSetDestroy(*cr, CRDestroy_Private));
33841b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
33941b6fd38SMatthew G. Knepley }
34041b6fd38SMatthew G. Knepley 
34197d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
3424b9ad928SBarry Smith {
343f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
344ce94432eSBarry Smith   MPI_Comm       comm;
3453aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
34610eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
34710167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
348f3fbd535SBarry Smith   PetscInt       i;
349f3fbd535SBarry Smith   PetscMPIInt    size;
350f3fbd535SBarry Smith   const char     *prefix;
351f3fbd535SBarry Smith   PC             ipc;
3523aeaf226SBarry Smith   PetscInt       n;
3534b9ad928SBarry Smith 
3544b9ad928SBarry Smith   PetscFunctionBegin;
3550700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
356548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
357548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
3589566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
3593aeaf226SBarry Smith   if (mglevels) {
36010eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
3613aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
3629566063dSJacob Faibussowitsch     PetscCall(PCReset_MG(pc));
3633aeaf226SBarry Smith     n    = mglevels[0]->levels;
3643aeaf226SBarry Smith     for (i=0; i<n; i++) {
3653aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
3669566063dSJacob Faibussowitsch         PetscCall(KSPDestroy(&mglevels[i]->smoothd));
3673aeaf226SBarry Smith       }
3689566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&mglevels[i]->smoothu));
3699566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&mglevels[i]->cr));
3709566063dSJacob Faibussowitsch       PetscCall(PetscFree(mglevels[i]));
3713aeaf226SBarry Smith     }
3729566063dSJacob Faibussowitsch     PetscCall(PetscFree(mg->levels));
3733aeaf226SBarry Smith   }
374f3fbd535SBarry Smith 
375f3fbd535SBarry Smith   mg->nlevels = levels;
376f3fbd535SBarry Smith 
3779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(levels,&mglevels));
3789566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*))));
379f3fbd535SBarry Smith 
3809566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc,&prefix));
381f3fbd535SBarry Smith 
382a9db87a2SMatthew G Knepley   mg->stageApply = 0;
383f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
3849566063dSJacob Faibussowitsch     PetscCall(PetscNewLog(pc,&mglevels[i]));
3852fa5cd67SKarl Rupp 
386f3fbd535SBarry Smith     mglevels[i]->level               = i;
387f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
38810eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
389336babb1SJed Brown     mg->default_smoothu              = 2;
390336babb1SJed Brown     mg->default_smoothd              = 2;
39163e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
39263e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
39363e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
39463e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
395f3fbd535SBarry Smith 
396f3fbd535SBarry Smith     if (comms) comm = comms[i];
397d5a8f7e6SBarry Smith     if (comm != MPI_COMM_NULL) {
3989566063dSJacob Faibussowitsch       PetscCall(KSPCreate(comm,&mglevels[i]->smoothd));
3999566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure));
4009566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i));
4019566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix));
4029566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level));
4030ee9a668SBarry Smith       if (i || levels == 1) {
4040ee9a668SBarry Smith         char tprefix[128];
4050ee9a668SBarry Smith 
4069566063dSJacob Faibussowitsch         PetscCall(KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV));
4079566063dSJacob Faibussowitsch         PetscCall(KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL));
4089566063dSJacob Faibussowitsch         PetscCall(KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE));
4099566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(mglevels[i]->smoothd,&ipc));
4109566063dSJacob Faibussowitsch         PetscCall(PCSetType(ipc,PCSOR));
4119566063dSJacob Faibussowitsch         PetscCall(KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd));
412f3fbd535SBarry Smith 
4139566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i));
4149566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix));
4150ee9a668SBarry Smith       } else {
4169566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_"));
417f3fbd535SBarry Smith 
4189dbfc187SHong Zhang         /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
4199566063dSJacob Faibussowitsch         PetscCall(KSPSetType(mglevels[0]->smoothd,KSPPREONLY));
4209566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(mglevels[0]->smoothd,&ipc));
4219566063dSJacob Faibussowitsch         PetscCallMPI(MPI_Comm_size(comm,&size));
422f3fbd535SBarry Smith         if (size > 1) {
4239566063dSJacob Faibussowitsch           PetscCall(PCSetType(ipc,PCREDUNDANT));
424f3fbd535SBarry Smith         } else {
4259566063dSJacob Faibussowitsch           PetscCall(PCSetType(ipc,PCLU));
426f3fbd535SBarry Smith         }
4279566063dSJacob Faibussowitsch         PetscCall(PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS));
428f3fbd535SBarry Smith       }
4299566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd));
430d5a8f7e6SBarry Smith     }
431f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
43231567311SBarry Smith     mg->rtol             = 0.0;
43331567311SBarry Smith     mg->abstol           = 0.0;
43431567311SBarry Smith     mg->dtol             = 0.0;
43531567311SBarry Smith     mg->ttol             = 0.0;
43631567311SBarry Smith     mg->cyclesperpcapply = 1;
437f3fbd535SBarry Smith   }
438f3fbd535SBarry Smith   mg->levels = mglevels;
4399566063dSJacob Faibussowitsch   PetscCall(PCMGSetType(pc,mgtype));
4404b9ad928SBarry Smith   PetscFunctionReturn(0);
4414b9ad928SBarry Smith }
4424b9ad928SBarry Smith 
44397d33e41SMatthew G. Knepley /*@C
44497d33e41SMatthew G. Knepley    PCMGSetLevels - Sets the number of levels to use with MG.
44597d33e41SMatthew G. Knepley    Must be called before any other MG routine.
44697d33e41SMatthew G. Knepley 
44797d33e41SMatthew G. Knepley    Logically Collective on PC
44897d33e41SMatthew G. Knepley 
44997d33e41SMatthew G. Knepley    Input Parameters:
45097d33e41SMatthew G. Knepley +  pc - the preconditioner context
45197d33e41SMatthew G. Knepley .  levels - the number of levels
45297d33e41SMatthew G. Knepley -  comms - optional communicators for each level; this is to allow solving the coarser problems
453d5a8f7e6SBarry Smith            on smaller sets of processes. For processes that are not included in the computation
45405552d4bSJunchao Zhang            you must pass MPI_COMM_NULL. Use comms = NULL to specify that all processes
45505552d4bSJunchao Zhang            should participate in each level of problem.
45697d33e41SMatthew G. Knepley 
45797d33e41SMatthew G. Knepley    Level: intermediate
45897d33e41SMatthew G. Knepley 
45997d33e41SMatthew G. Knepley    Notes:
46097d33e41SMatthew G. Knepley      If the number of levels is one then the multigrid uses the -mg_levels prefix
46197d33e41SMatthew G. Knepley      for setting the level options rather than the -mg_coarse prefix.
46297d33e41SMatthew G. Knepley 
463d5a8f7e6SBarry Smith      You can free the information in comms after this routine is called.
464d5a8f7e6SBarry Smith 
465d5a8f7e6SBarry Smith      The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level
466d5a8f7e6SBarry Smith      are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on
467d5a8f7e6SBarry Smith      the two levels, rank 0 in the original communicator will pass in an array of 2 communicators
468d5a8f7e6SBarry Smith      of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators
469d5a8f7e6SBarry Smith      the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate
470d5a8f7e6SBarry Smith      in the coarse grid solve.
471d5a8f7e6SBarry Smith 
472d5a8f7e6SBarry Smith      Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one
473d5a8f7e6SBarry Smith      must take special care in providing the restriction and interpolation operation. We recommend
474d5a8f7e6SBarry Smith      providing these as two step operations; first perform a standard restriction or interpolation on
475d5a8f7e6SBarry Smith      the full number of ranks for that level and then use an MPI call to copy the resulting vector
47605552d4bSJunchao Zhang      array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, note in both
477d5a8f7e6SBarry Smith      cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and
478d5a8f7e6SBarry Smith      recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries.
479d5a8f7e6SBarry Smith 
48005552d4bSJunchao Zhang    Fortran Notes:
48105552d4bSJunchao Zhang      Use comms = PETSC_NULL_MPI_COMM as the equivalent of NULL in the C interface. Note PETSC_NULL_MPI_COMM
48205552d4bSJunchao Zhang      is not MPI_COMM_NULL. It is more like PETSC_NULL_INTEGER, PETSC_NULL_REAL etc.
483d5a8f7e6SBarry Smith 
48497d33e41SMatthew G. Knepley .seealso: PCMGSetType(), PCMGGetLevels()
48597d33e41SMatthew G. Knepley @*/
48697d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
48797d33e41SMatthew G. Knepley {
48897d33e41SMatthew G. Knepley   PetscFunctionBegin;
48997d33e41SMatthew G. Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
49097d33e41SMatthew G. Knepley   if (comms) PetscValidPointer(comms,3);
4919566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms)));
49297d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
49397d33e41SMatthew G. Knepley }
49497d33e41SMatthew G. Knepley 
495c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
496c07bf074SBarry Smith {
497a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
498a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
499a06653b4SBarry Smith   PetscInt       i,n;
500c07bf074SBarry Smith 
501c07bf074SBarry Smith   PetscFunctionBegin;
5029566063dSJacob Faibussowitsch   PetscCall(PCReset_MG(pc));
503a06653b4SBarry Smith   if (mglevels) {
504a06653b4SBarry Smith     n = mglevels[0]->levels;
505a06653b4SBarry Smith     for (i=0; i<n; i++) {
506a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
5079566063dSJacob Faibussowitsch         PetscCall(KSPDestroy(&mglevels[i]->smoothd));
508a06653b4SBarry Smith       }
5099566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&mglevels[i]->smoothu));
5109566063dSJacob Faibussowitsch       PetscCall(KSPDestroy(&mglevels[i]->cr));
5119566063dSJacob Faibussowitsch       PetscCall(PetscFree(mglevels[i]));
512a06653b4SBarry Smith     }
5139566063dSJacob Faibussowitsch     PetscCall(PetscFree(mg->levels));
514a06653b4SBarry Smith   }
5159566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
5169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL));
5179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL));
518f3fbd535SBarry Smith   PetscFunctionReturn(0);
519f3fbd535SBarry Smith }
520f3fbd535SBarry Smith 
521f3fbd535SBarry Smith /*
522f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
523f3fbd535SBarry Smith              or full cycle of multigrid.
524f3fbd535SBarry Smith 
525f3fbd535SBarry Smith   Note:
526f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
527f3fbd535SBarry Smith */
52830b0564aSStefano Zampini static PetscErrorCode PCApply_MG_Internal(PC pc,Vec b,Vec x,Mat B,Mat X,PetscBool transpose)
529f3fbd535SBarry Smith {
530f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
531f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
532391689abSStefano Zampini   PC             tpc;
533f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
53430b0564aSStefano Zampini   PetscBool      changeu,changed,matapp;
535f3fbd535SBarry Smith 
536f3fbd535SBarry Smith   PetscFunctionBegin;
53730b0564aSStefano Zampini   matapp = (PetscBool)(B && X);
5389566063dSJacob Faibussowitsch   if (mg->stageApply) PetscCall(PetscLogStagePush(mg->stageApply));
539e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
540298cc208SBarry Smith   for (i=0; i<levels; i++) {
541e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
5429566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL));
5439566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A));
544e1d8e5deSBarry Smith     }
545e1d8e5deSBarry Smith   }
546e1d8e5deSBarry Smith 
5479566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[levels-1]->smoothd,&tpc));
5489566063dSJacob Faibussowitsch   PetscCall(PCPreSolveChangeRHS(tpc,&changed));
5499566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[levels-1]->smoothu,&tpc));
5509566063dSJacob Faibussowitsch   PetscCall(PCPreSolveChangeRHS(tpc,&changeu));
551391689abSStefano Zampini   if (!changeu && !changed) {
55230b0564aSStefano Zampini     if (matapp) {
5539566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[levels-1]->B));
55430b0564aSStefano Zampini       mglevels[levels-1]->B = B;
55530b0564aSStefano Zampini     } else {
5569566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&mglevels[levels-1]->b));
557f3fbd535SBarry Smith       mglevels[levels-1]->b = b;
55830b0564aSStefano Zampini     }
559391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
56030b0564aSStefano Zampini     if (matapp) {
56130b0564aSStefano Zampini       if (mglevels[levels-1]->B) {
56230b0564aSStefano Zampini         PetscInt  N1,N2;
56330b0564aSStefano Zampini         PetscBool flg;
56430b0564aSStefano Zampini 
5659566063dSJacob Faibussowitsch         PetscCall(MatGetSize(mglevels[levels-1]->B,NULL,&N1));
5669566063dSJacob Faibussowitsch         PetscCall(MatGetSize(B,NULL,&N2));
5679566063dSJacob Faibussowitsch         PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels-1]->B,((PetscObject)B)->type_name,&flg));
56830b0564aSStefano Zampini         if (N1 != N2 || !flg) {
5699566063dSJacob Faibussowitsch           PetscCall(MatDestroy(&mglevels[levels-1]->B));
57030b0564aSStefano Zampini         }
57130b0564aSStefano Zampini       }
57230b0564aSStefano Zampini       if (!mglevels[levels-1]->B) {
5739566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&mglevels[levels-1]->B));
57430b0564aSStefano Zampini       } else {
5759566063dSJacob Faibussowitsch         PetscCall(MatCopy(B,mglevels[levels-1]->B,SAME_NONZERO_PATTERN));
57630b0564aSStefano Zampini       }
57730b0564aSStefano Zampini     } else {
578391689abSStefano Zampini       if (!mglevels[levels-1]->b) {
579391689abSStefano Zampini         Vec *vec;
580391689abSStefano Zampini 
5819566063dSJacob Faibussowitsch         PetscCall(KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL));
582391689abSStefano Zampini         mglevels[levels-1]->b = *vec;
5839566063dSJacob Faibussowitsch         PetscCall(PetscFree(vec));
584391689abSStefano Zampini       }
5859566063dSJacob Faibussowitsch       PetscCall(VecCopy(b,mglevels[levels-1]->b));
586391689abSStefano Zampini     }
58730b0564aSStefano Zampini   }
58830b0564aSStefano Zampini   if (matapp) { mglevels[levels-1]->X = X; }
58930b0564aSStefano Zampini   else { mglevels[levels-1]->x = x; }
59030b0564aSStefano Zampini 
59130b0564aSStefano Zampini   /* If coarser Xs are present, it means we have already block applied the PC at least once
59230b0564aSStefano Zampini      Reset operators if sizes/type do no match */
59330b0564aSStefano Zampini   if (matapp && levels > 1 && mglevels[levels-2]->X) {
59430b0564aSStefano Zampini     PetscInt  Xc,Bc;
59530b0564aSStefano Zampini     PetscBool flg;
59630b0564aSStefano Zampini 
5979566063dSJacob Faibussowitsch     PetscCall(MatGetSize(mglevels[levels-2]->X,NULL,&Xc));
5989566063dSJacob Faibussowitsch     PetscCall(MatGetSize(mglevels[levels-1]->B,NULL,&Bc));
5999566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels-2]->X,((PetscObject)mglevels[levels-1]->X)->type_name,&flg));
60030b0564aSStefano Zampini     if (Xc != Bc || !flg) {
6019566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&mglevels[levels-1]->R));
60230b0564aSStefano Zampini       for (i=0;i<levels-1;i++) {
6039566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&mglevels[i]->R));
6049566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&mglevels[i]->B));
6059566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&mglevels[i]->X));
60630b0564aSStefano Zampini       }
60730b0564aSStefano Zampini     }
60830b0564aSStefano Zampini   }
609391689abSStefano Zampini 
61031567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
6119566063dSJacob Faibussowitsch     if (matapp) PetscCall(MatZeroEntries(X));
6129566063dSJacob Faibussowitsch     else PetscCall(VecZeroEntries(x));
61331567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
6149566063dSJacob Faibussowitsch       PetscCall(PCMGMCycle_Private(pc,mglevels+levels-1,transpose,matapp,NULL));
615f3fbd535SBarry Smith     }
6162fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
6179566063dSJacob Faibussowitsch     PetscCall(PCMGACycle_Private(pc,mglevels,transpose,matapp));
6182fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
6199566063dSJacob Faibussowitsch     PetscCall(PCMGKCycle_Private(pc,mglevels,transpose,matapp));
6202fa5cd67SKarl Rupp   } else {
6219566063dSJacob Faibussowitsch     PetscCall(PCMGFCycle_Private(pc,mglevels,transpose,matapp));
622f3fbd535SBarry Smith   }
6239566063dSJacob Faibussowitsch   if (mg->stageApply) PetscCall(PetscLogStagePop());
62430b0564aSStefano Zampini   if (!changeu && !changed) {
62530b0564aSStefano Zampini     if (matapp) { mglevels[levels-1]->B = NULL; }
62630b0564aSStefano Zampini     else { mglevels[levels-1]->b = NULL; }
62730b0564aSStefano Zampini   }
628f3fbd535SBarry Smith   PetscFunctionReturn(0);
629f3fbd535SBarry Smith }
630f3fbd535SBarry Smith 
631fcb023d4SJed Brown static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
632fcb023d4SJed Brown {
633fcb023d4SJed Brown   PetscFunctionBegin;
6349566063dSJacob Faibussowitsch   PetscCall(PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_FALSE));
635fcb023d4SJed Brown   PetscFunctionReturn(0);
636fcb023d4SJed Brown }
637fcb023d4SJed Brown 
638fcb023d4SJed Brown static PetscErrorCode PCApplyTranspose_MG(PC pc,Vec b,Vec x)
639fcb023d4SJed Brown {
640fcb023d4SJed Brown   PetscFunctionBegin;
6419566063dSJacob Faibussowitsch   PetscCall(PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_TRUE));
64230b0564aSStefano Zampini   PetscFunctionReturn(0);
64330b0564aSStefano Zampini }
64430b0564aSStefano Zampini 
64530b0564aSStefano Zampini static PetscErrorCode PCMatApply_MG(PC pc,Mat b,Mat x)
64630b0564aSStefano Zampini {
64730b0564aSStefano Zampini   PetscFunctionBegin;
6489566063dSJacob Faibussowitsch   PetscCall(PCApply_MG_Internal(pc,NULL,NULL,b,x,PETSC_FALSE));
649fcb023d4SJed Brown   PetscFunctionReturn(0);
650fcb023d4SJed Brown }
651f3fbd535SBarry Smith 
6524416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
653f3fbd535SBarry Smith {
654710315b6SLawrence Mitchell   PetscInt         levels,cycles;
655f3b08a26SMatthew G. Knepley   PetscBool        flg, flg2;
656f3fbd535SBarry Smith   PC_MG            *mg = (PC_MG*)pc->data;
6573d3eaba7SBarry Smith   PC_MG_Levels     **mglevels;
658f3fbd535SBarry Smith   PCMGType         mgtype;
659f3fbd535SBarry Smith   PCMGCycleType    mgctype;
6602134b1e4SBarry Smith   PCMGGalerkinType gtype;
661f3fbd535SBarry Smith 
662f3fbd535SBarry Smith   PetscFunctionBegin;
6631d743356SStefano Zampini   levels = PetscMax(mg->nlevels,1);
6649566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"Multigrid options"));
6659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg));
6661a97d7b6SStefano Zampini   if (!flg && !mg->levels && pc->dm) {
6679566063dSJacob Faibussowitsch     PetscCall(DMGetRefineLevel(pc->dm,&levels));
668eb3f98d2SBarry Smith     levels++;
6693aeaf226SBarry Smith     mg->usedmfornumberoflevels = PETSC_TRUE;
670eb3f98d2SBarry Smith   }
6719566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels(pc,levels,NULL));
6723aeaf226SBarry Smith   mglevels = mg->levels;
6733aeaf226SBarry Smith 
674f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
6759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg));
676f3fbd535SBarry Smith   if (flg) {
6779566063dSJacob Faibussowitsch     PetscCall(PCMGSetCycleType(pc,mgctype));
6782fa5cd67SKarl Rupp   }
6792134b1e4SBarry Smith   gtype = mg->galerkin;
6809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg));
6812134b1e4SBarry Smith   if (flg) {
6829566063dSJacob Faibussowitsch     PetscCall(PCMGSetGalerkin(pc,gtype));
683f3fbd535SBarry Smith   }
684f3b08a26SMatthew G. Knepley   flg2 = PETSC_FALSE;
6859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg));
6869566063dSJacob Faibussowitsch   if (flg) PetscCall(PCMGSetAdaptInterpolation(pc, flg2));
6879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg));
6889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_mg_adapt_interp_coarse_space","Type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector","PCMGSetAdaptCoarseSpaceType",PCMGCoarseSpaceTypes,(PetscEnum)mg->coarseSpaceType,(PetscEnum*)&mg->coarseSpaceType,&flg));
6899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg));
69041b6fd38SMatthew G. Knepley   flg2 = PETSC_FALSE;
6919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg));
6929566063dSJacob Faibussowitsch   if (flg) PetscCall(PCMGSetAdaptCR(pc, flg2));
693f56b1265SBarry Smith   flg = PETSC_FALSE;
6949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL));
695f442ab6aSBarry Smith   if (flg) {
6969566063dSJacob Faibussowitsch     PetscCall(PCMGSetDistinctSmoothUp(pc));
697f442ab6aSBarry Smith   }
69831567311SBarry Smith   mgtype = mg->am;
6999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg));
700f3fbd535SBarry Smith   if (flg) {
7019566063dSJacob Faibussowitsch     PetscCall(PCMGSetType(pc,mgtype));
702f3fbd535SBarry Smith   }
70331567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
7049566063dSJacob Faibussowitsch     PetscCall(PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg));
705f3fbd535SBarry Smith     if (flg) {
7069566063dSJacob Faibussowitsch       PetscCall(PCMGMultiplicativeSetCycles(pc,cycles));
707f3fbd535SBarry Smith     }
708f3fbd535SBarry Smith   }
709f3fbd535SBarry Smith   flg  = PETSC_FALSE;
7109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL));
711f3fbd535SBarry Smith   if (flg) {
712f3fbd535SBarry Smith     PetscInt i;
713f3fbd535SBarry Smith     char     eventname[128];
7141a97d7b6SStefano Zampini 
715f3fbd535SBarry Smith     levels = mglevels[0]->levels;
716f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
717f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
7189566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup));
719f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
7209566063dSJacob Faibussowitsch       PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve));
721f3fbd535SBarry Smith       if (i) {
722f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
7239566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual));
724f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
7259566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict));
726f3fbd535SBarry Smith       }
727f3fbd535SBarry Smith     }
728eec35531SMatthew G Knepley 
729ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
730eec35531SMatthew G Knepley     {
731eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
732eec35531SMatthew G Knepley       PetscStageLog stageLog;
733eec35531SMatthew G Knepley       PetscInt      st;
734eec35531SMatthew G Knepley 
7359566063dSJacob Faibussowitsch       PetscCall(PetscLogGetStageLog(&stageLog));
736eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
737eec35531SMatthew G Knepley         PetscBool same;
738eec35531SMatthew G Knepley 
7399566063dSJacob Faibussowitsch         PetscCall(PetscStrcmp(stageLog->stageInfo[st].name, sname, &same));
7402fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
741eec35531SMatthew G Knepley       }
742eec35531SMatthew G Knepley       if (!mg->stageApply) {
7439566063dSJacob Faibussowitsch         PetscCall(PetscLogStageRegister(sname, &mg->stageApply));
744eec35531SMatthew G Knepley       }
745eec35531SMatthew G Knepley     }
746ec60ca73SSatish Balay #endif
747f3fbd535SBarry Smith   }
7489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
749f3b08a26SMatthew G. Knepley   /* Check option consistency */
7509566063dSJacob Faibussowitsch   PetscCall(PCMGGetGalerkin(pc, &gtype));
7519566063dSJacob Faibussowitsch   PetscCall(PCMGGetAdaptInterpolation(pc, &flg));
7522c71b3e2SJacob Faibussowitsch   PetscCheckFalse(flg && (gtype >= PC_MG_GALERKIN_NONE),PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
753f3fbd535SBarry Smith   PetscFunctionReturn(0);
754f3fbd535SBarry Smith }
755f3fbd535SBarry Smith 
7560a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL};
7570a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL};
7580a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL};
759f3b08a26SMatthew G. Knepley const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL};
760f3fbd535SBarry Smith 
7619804daf3SBarry Smith #include <petscdraw.h>
762f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
763f3fbd535SBarry Smith {
764f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
765f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
766e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
767219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
768f3fbd535SBarry Smith 
769f3fbd535SBarry Smith   PetscFunctionBegin;
7709566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
7719566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary));
7729566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw));
773f3fbd535SBarry Smith   if (iascii) {
774e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
7759566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename));
77631567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
7779566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply));
778f3fbd535SBarry Smith     }
7792134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
7809566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n"));
7812134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
7829566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n"));
7832134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
7849566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n"));
7852134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
7869566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n"));
7874f66f45eSBarry Smith     } else {
7889566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n"));
789f3fbd535SBarry Smith     }
7905adeb434SBarry Smith     if (mg->view) {
7919566063dSJacob Faibussowitsch       PetscCall((*mg->view)(pc,viewer));
7925adeb434SBarry Smith     }
793f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
794f3fbd535SBarry Smith       if (!i) {
7959566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i));
796f3fbd535SBarry Smith       } else {
7979566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i));
798f3fbd535SBarry Smith       }
7999566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
8009566063dSJacob Faibussowitsch       PetscCall(KSPView(mglevels[i]->smoothd,viewer));
8019566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
802f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
8039566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n"));
804f3fbd535SBarry Smith       } else if (i) {
8059566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i));
8069566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
8079566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothu,viewer));
8089566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
809f3fbd535SBarry Smith       }
81041b6fd38SMatthew G. Knepley       if (i && mglevels[i]->cr) {
8119566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer,"CR solver on level %D -------------------------------\n",i));
8129566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushTab(viewer));
8139566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->cr,viewer));
8149566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopTab(viewer));
81541b6fd38SMatthew G. Knepley       }
816f3fbd535SBarry Smith     }
8175b0b0462SBarry Smith   } else if (isbinary) {
8185b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
8199566063dSJacob Faibussowitsch       PetscCall(KSPView(mglevels[i]->smoothd,viewer));
8205b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
8219566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothu,viewer));
8225b0b0462SBarry Smith       }
8235b0b0462SBarry Smith     }
824219b1045SBarry Smith   } else if (isdraw) {
825219b1045SBarry Smith     PetscDraw draw;
826b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
8279566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw));
8289566063dSJacob Faibussowitsch     PetscCall(PetscDrawGetCurrentPoint(draw,&x,&y));
8299566063dSJacob Faibussowitsch     PetscCall(PetscDrawStringGetSize(draw,NULL,&th));
830219b1045SBarry Smith     bottom = y - th;
831219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
832b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
8339566063dSJacob Faibussowitsch         PetscCall(PetscDrawPushCurrentPoint(draw,x,bottom));
8349566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothd,viewer));
8359566063dSJacob Faibussowitsch         PetscCall(PetscDrawPopCurrentPoint(draw));
836b4375e8dSPeter Brune       } else {
837b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
8389566063dSJacob Faibussowitsch         PetscCall(PetscDrawPushCurrentPoint(draw,x+w,bottom));
8399566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothd,viewer));
8409566063dSJacob Faibussowitsch         PetscCall(PetscDrawPopCurrentPoint(draw));
8419566063dSJacob Faibussowitsch         PetscCall(PetscDrawPushCurrentPoint(draw,x-w,bottom));
8429566063dSJacob Faibussowitsch         PetscCall(KSPView(mglevels[i]->smoothu,viewer));
8439566063dSJacob Faibussowitsch         PetscCall(PetscDrawPopCurrentPoint(draw));
844b4375e8dSPeter Brune       }
8459566063dSJacob Faibussowitsch       PetscCall(PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL));
8461cd381d6SBarry Smith       bottom -= th;
847219b1045SBarry Smith     }
8485b0b0462SBarry Smith   }
849f3fbd535SBarry Smith   PetscFunctionReturn(0);
850f3fbd535SBarry Smith }
851f3fbd535SBarry Smith 
852af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
853cab2e9ccSBarry Smith 
854f3fbd535SBarry Smith /*
855f3fbd535SBarry Smith     Calls setup for the KSP on each level
856f3fbd535SBarry Smith */
857f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
858f3fbd535SBarry Smith {
859f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
860f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
8611a97d7b6SStefano Zampini   PetscInt       i,n;
86298e04984SBarry Smith   PC             cpc;
8633db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
864f3fbd535SBarry Smith   Mat            dA,dB;
865f3fbd535SBarry Smith   Vec            tvec;
866218a07d4SBarry Smith   DM             *dms;
8670a545947SLisandro Dalcin   PetscViewer    viewer = NULL;
86841b6fd38SMatthew G. Knepley   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE;
869f3fbd535SBarry Smith 
870f3fbd535SBarry Smith   PetscFunctionBegin;
87128b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
8721a97d7b6SStefano Zampini   n = mglevels[0]->levels;
87301bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
8743aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
8753aeaf226SBarry Smith     PetscInt levels;
8769566063dSJacob Faibussowitsch     PetscCall(DMGetRefineLevel(pc->dm,&levels));
8773aeaf226SBarry Smith     levels++;
8783aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
8799566063dSJacob Faibussowitsch       PetscCall(PCMGSetLevels(pc,levels,NULL));
8803aeaf226SBarry Smith       n        = levels;
8819566063dSJacob Faibussowitsch       PetscCall(PCSetFromOptions(pc)); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
8823aeaf226SBarry Smith       mglevels = mg->levels;
8833aeaf226SBarry Smith     }
8843aeaf226SBarry Smith   }
8859566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(mglevels[0]->smoothd,&cpc));
8863aeaf226SBarry Smith 
887f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
888f3fbd535SBarry Smith   /* so use those from global PC */
889f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
8909566063dSJacob Faibussowitsch   PetscCall(KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset));
89198e04984SBarry Smith   if (opsset) {
89298e04984SBarry Smith     Mat mmat;
8939566063dSJacob Faibussowitsch     PetscCall(KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat));
89498e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
89598e04984SBarry Smith   }
8964a5f07a5SMark F. Adams 
89741b6fd38SMatthew G. Knepley   /* Create CR solvers */
8989566063dSJacob Faibussowitsch   PetscCall(PCMGGetAdaptCR(pc, &doCR));
89941b6fd38SMatthew G. Knepley   if (doCR) {
90041b6fd38SMatthew G. Knepley     const char *prefix;
90141b6fd38SMatthew G. Knepley 
9029566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
90341b6fd38SMatthew G. Knepley     for (i = 1; i < n; ++i) {
90441b6fd38SMatthew G. Knepley       PC   ipc, cr;
90541b6fd38SMatthew G. Knepley       char crprefix[128];
90641b6fd38SMatthew G. Knepley 
9079566063dSJacob Faibussowitsch       PetscCall(KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr));
9089566063dSJacob Faibussowitsch       PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE));
9099566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i));
9109566063dSJacob Faibussowitsch       PetscCall(KSPSetOptionsPrefix(mglevels[i]->cr, prefix));
9119566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level));
9129566063dSJacob Faibussowitsch       PetscCall(KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV));
9139566063dSJacob Faibussowitsch       PetscCall(KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL));
9149566063dSJacob Faibussowitsch       PetscCall(KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED));
9159566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(mglevels[i]->cr, &ipc));
91641b6fd38SMatthew G. Knepley 
9179566063dSJacob Faibussowitsch       PetscCall(PCSetType(ipc, PCCOMPOSITE));
9189566063dSJacob Faibussowitsch       PetscCall(PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE));
9199566063dSJacob Faibussowitsch       PetscCall(PCCompositeAddPCType(ipc, PCSOR));
9209566063dSJacob Faibussowitsch       PetscCall(CreateCR_Private(pc, i, &cr));
9219566063dSJacob Faibussowitsch       PetscCall(PCCompositeAddPC(ipc, cr));
9229566063dSJacob Faibussowitsch       PetscCall(PCDestroy(&cr));
92341b6fd38SMatthew G. Knepley 
9249566063dSJacob Faibussowitsch       PetscCall(KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd));
9259566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE));
9269566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i));
9279566063dSJacob Faibussowitsch       PetscCall(KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix));
9289566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr));
92941b6fd38SMatthew G. Knepley     }
93041b6fd38SMatthew G. Knepley   }
93141b6fd38SMatthew G. Knepley 
9324a5f07a5SMark F. Adams   if (!opsset) {
9339566063dSJacob Faibussowitsch     PetscCall(PCGetUseAmat(pc,&use_amat));
9344a5f07a5SMark F. Adams     if (use_amat) {
9359566063dSJacob 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"));
9369566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat));
93769858f1bSStefano Zampini     } else {
9389566063dSJacob 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"));
9399566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat));
9404a5f07a5SMark F. Adams     }
9414a5f07a5SMark F. Adams   }
942f3fbd535SBarry Smith 
9439c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
9449c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
9459c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
9469c7ffb39SBarry Smith       continue;
9479c7ffb39SBarry Smith     }
9489c7ffb39SBarry Smith   }
9492134b1e4SBarry Smith 
9509566063dSJacob Faibussowitsch   PetscCall(KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB));
9512134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
9522134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
9532134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
9542134b1e4SBarry Smith   }
9552134b1e4SBarry Smith 
9569c7ffb39SBarry Smith   /*
9579c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
9589c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
9599c7ffb39SBarry Smith   */
9602134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
9612d2b81a6SBarry Smith         /* construct the interpolation from the DMs */
9622e499ae9SBarry Smith     Mat p;
96373250ac0SBarry Smith     Vec rscale;
9649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n,&dms));
965218a07d4SBarry Smith     dms[n-1] = pc->dm;
966ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
9679566063dSJacob Faibussowitsch     for (i=n-2; i>-1; i--) PetscCall(DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]));
96841fce8fdSHong Zhang         /*
96941fce8fdSHong Zhang            Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
97041fce8fdSHong Zhang            Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
97141fce8fdSHong Zhang            But it is safe to use -dm_mat_type.
97241fce8fdSHong Zhang 
97341fce8fdSHong Zhang            The mat type should not be hardcoded like this, we need to find a better way.
9749566063dSJacob Faibussowitsch     PetscCall(DMSetMatType(dms[0],MATAIJ));
97541fce8fdSHong Zhang     */
976218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
977942e3340SBarry Smith       DMKSP     kdm;
978eab52d2bSLawrence Mitchell       PetscBool dmhasrestrict, dmhasinject;
9799566063dSJacob Faibussowitsch       PetscCall(KSPSetDM(mglevels[i]->smoothd,dms[i]));
9809566063dSJacob Faibussowitsch       if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE));
981c27ee7a3SPatrick Farrell       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
9829566063dSJacob Faibussowitsch         PetscCall(KSPSetDM(mglevels[i]->smoothu,dms[i]));
9839566063dSJacob Faibussowitsch         if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE));
984c27ee7a3SPatrick Farrell       }
98541b6fd38SMatthew G. Knepley       if (mglevels[i]->cr) {
9869566063dSJacob Faibussowitsch         PetscCall(KSPSetDM(mglevels[i]->cr,dms[i]));
9879566063dSJacob Faibussowitsch         if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE));
98841b6fd38SMatthew G. Knepley       }
9899566063dSJacob Faibussowitsch       PetscCall(DMGetDMKSPWrite(dms[i],&kdm));
990d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
991d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
9920298fd71SBarry Smith       kdm->ops->computerhs = NULL;
9930298fd71SBarry Smith       kdm->rhsctx          = NULL;
99424c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
9959566063dSJacob Faibussowitsch         PetscCall(DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale));
9969566063dSJacob Faibussowitsch         PetscCall(PCMGSetInterpolation(pc,i+1,p));
9979566063dSJacob Faibussowitsch         if (rscale) PetscCall(PCMGSetRScale(pc,i+1,rscale));
9989566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&rscale));
9999566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&p));
1000218a07d4SBarry Smith       }
10019566063dSJacob Faibussowitsch       PetscCall(DMHasCreateRestriction(dms[i],&dmhasrestrict));
10023ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct) {
10039566063dSJacob Faibussowitsch         PetscCall(DMCreateRestriction(dms[i],dms[i+1],&p));
10049566063dSJacob Faibussowitsch         PetscCall(PCMGSetRestriction(pc,i+1,p));
10059566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&p));
10063ad4599aSBarry Smith       }
10079566063dSJacob Faibussowitsch       PetscCall(DMHasCreateInjection(dms[i],&dmhasinject));
1008eab52d2bSLawrence Mitchell       if (dmhasinject && !mglevels[i+1]->inject) {
10099566063dSJacob Faibussowitsch         PetscCall(DMCreateInjection(dms[i],dms[i+1],&p));
10109566063dSJacob Faibussowitsch         PetscCall(PCMGSetInjection(pc,i+1,p));
10119566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&p));
1012eab52d2bSLawrence Mitchell       }
101324c3aa18SBarry Smith     }
10142d2b81a6SBarry Smith 
10159566063dSJacob Faibussowitsch     for (i=n-2; i>-1; i--) PetscCall(DMDestroy(&dms[i]));
10169566063dSJacob Faibussowitsch     PetscCall(PetscFree(dms));
1017ce4cda84SJed Brown   }
1018cab2e9ccSBarry Smith 
1019ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
10202134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
10219566063dSJacob Faibussowitsch     PetscCall(KSPSetDM(mglevels[n-1]->smoothd,pc->dm));
10229566063dSJacob Faibussowitsch     PetscCall(KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE));
1023c27ee7a3SPatrick Farrell     if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) {
10249566063dSJacob Faibussowitsch       PetscCall(KSPSetDM(mglevels[n-1]->smoothu,pc->dm));
10259566063dSJacob Faibussowitsch       PetscCall(KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE));
1026c27ee7a3SPatrick Farrell     }
102741b6fd38SMatthew G. Knepley     if (mglevels[n-1]->cr) {
10289566063dSJacob Faibussowitsch       PetscCall(KSPSetDM(mglevels[n-1]->cr,pc->dm));
10299566063dSJacob Faibussowitsch       PetscCall(KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE));
103041b6fd38SMatthew G. Knepley     }
1031218a07d4SBarry Smith   }
1032218a07d4SBarry Smith 
10332134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
10342134b1e4SBarry Smith     Mat       A,B;
10352134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
10362134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
10372134b1e4SBarry Smith 
10382134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
10392134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
10402134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
1041f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
1042*7827d75bSBarry 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");
1043b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
10449566063dSJacob Faibussowitsch         PetscCall(PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct));
1045b9367914SBarry Smith       }
1046b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
10479566063dSJacob Faibussowitsch         PetscCall(PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate));
1048b9367914SBarry Smith       }
10492134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
10509566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd,&A,&B));
10512134b1e4SBarry Smith       }
10522134b1e4SBarry Smith       if (doA) {
10539566063dSJacob Faibussowitsch         PetscCall(MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A));
10542134b1e4SBarry Smith       }
10552134b1e4SBarry Smith       if (doB) {
10569566063dSJacob Faibussowitsch         PetscCall(MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B));
1057b9367914SBarry Smith       }
10582134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
10592134b1e4SBarry Smith       if (!doA && dAeqdB) {
10609566063dSJacob Faibussowitsch         if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)B));
10612134b1e4SBarry Smith         A = B;
10622134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
10639566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd,&A,NULL));
10649566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)A));
1065b9367914SBarry Smith       }
10662134b1e4SBarry Smith       if (!doB && dAeqdB) {
10679566063dSJacob Faibussowitsch         if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)A));
10682134b1e4SBarry Smith         B = A;
10692134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
10709566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd,NULL,&B));
10719566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)B));
10727d537d92SJed Brown       }
10732134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
10749566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[i]->smoothd,A,B));
10759566063dSJacob Faibussowitsch         PetscCall(PetscObjectDereference((PetscObject)A));
10769566063dSJacob Faibussowitsch         PetscCall(PetscObjectDereference((PetscObject)B));
10772134b1e4SBarry Smith       }
10782134b1e4SBarry Smith       dA = A;
1079f3fbd535SBarry Smith       dB = B;
1080f3fbd535SBarry Smith     }
1081f3fbd535SBarry Smith   }
1082f3b08a26SMatthew G. Knepley 
1083f3b08a26SMatthew G. Knepley   /* Adapt interpolation matrices */
1084f3b08a26SMatthew G. Knepley   if (mg->adaptInterpolation) {
1085f3b08a26SMatthew G. Knepley     mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */
1086f3b08a26SMatthew G. Knepley     for (i = 0; i < n; ++i) {
10879566063dSJacob Faibussowitsch       PetscCall(PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace));
10889566063dSJacob Faibussowitsch       if (i) PetscCall(PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mg->Nc, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace));
1089f3b08a26SMatthew G. Knepley     }
1090f3b08a26SMatthew G. Knepley     for (i = n-2; i > -1; --i) {
10919566063dSJacob Faibussowitsch       PetscCall(PCMGRecomputeLevelOperators_Internal(pc, i));
1092f3b08a26SMatthew G. Knepley     }
1093f3b08a26SMatthew G. Knepley   }
1094f3b08a26SMatthew G. Knepley 
10952134b1e4SBarry Smith   if (needRestricts && pc->dm) {
1096caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
1097ccceb50cSJed Brown       DM  dmfine,dmcoarse;
1098caa4e7f2SJed Brown       Mat Restrict,Inject;
1099caa4e7f2SJed Brown       Vec rscale;
11009566063dSJacob Faibussowitsch       PetscCall(KSPGetDM(mglevels[i+1]->smoothd,&dmfine));
11019566063dSJacob Faibussowitsch       PetscCall(KSPGetDM(mglevels[i]->smoothd,&dmcoarse));
11029566063dSJacob Faibussowitsch       PetscCall(PCMGGetRestriction(pc,i+1,&Restrict));
11039566063dSJacob Faibussowitsch       PetscCall(PCMGGetRScale(pc,i+1,&rscale));
11049566063dSJacob Faibussowitsch       PetscCall(PCMGGetInjection(pc,i+1,&Inject));
11059566063dSJacob Faibussowitsch       PetscCall(DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse));
1106caa4e7f2SJed Brown     }
1107caa4e7f2SJed Brown   }
1108f3fbd535SBarry Smith 
1109f3fbd535SBarry Smith   if (!pc->setupcalled) {
1110f3fbd535SBarry Smith     for (i=0; i<n; i++) {
11119566063dSJacob Faibussowitsch       PetscCall(KSPSetFromOptions(mglevels[i]->smoothd));
1112f3fbd535SBarry Smith     }
1113f3fbd535SBarry Smith     for (i=1; i<n; i++) {
1114f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
11159566063dSJacob Faibussowitsch         PetscCall(KSPSetFromOptions(mglevels[i]->smoothu));
1116f3fbd535SBarry Smith       }
111741b6fd38SMatthew G. Knepley       if (mglevels[i]->cr) {
11189566063dSJacob Faibussowitsch         PetscCall(KSPSetFromOptions(mglevels[i]->cr));
111941b6fd38SMatthew G. Knepley       }
1120f3fbd535SBarry Smith     }
11213ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
1122f3fbd535SBarry Smith     for (i=1; i<n; i++) {
11239566063dSJacob Faibussowitsch       PetscCall(PCMGGetInterpolation(pc,i,NULL));
11249566063dSJacob Faibussowitsch       PetscCall(PCMGGetRestriction(pc,i,NULL));
1125f3fbd535SBarry Smith     }
1126f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
1127f3fbd535SBarry Smith       if (!mglevels[i]->b) {
1128f3fbd535SBarry Smith         Vec *vec;
11299566063dSJacob Faibussowitsch         PetscCall(KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL));
11309566063dSJacob Faibussowitsch         PetscCall(PCMGSetRhs(pc,i,*vec));
11319566063dSJacob Faibussowitsch         PetscCall(VecDestroy(vec));
11329566063dSJacob Faibussowitsch         PetscCall(PetscFree(vec));
1133f3fbd535SBarry Smith       }
1134f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
11359566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b,&tvec));
11369566063dSJacob Faibussowitsch         PetscCall(PCMGSetR(pc,i,tvec));
11379566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&tvec));
1138f3fbd535SBarry Smith       }
1139f3fbd535SBarry Smith       if (!mglevels[i]->x) {
11409566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b,&tvec));
11419566063dSJacob Faibussowitsch         PetscCall(PCMGSetX(pc,i,tvec));
11429566063dSJacob Faibussowitsch         PetscCall(VecDestroy(&tvec));
1143f3fbd535SBarry Smith       }
114441b6fd38SMatthew G. Knepley       if (doCR) {
11459566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b,&mglevels[i]->crx));
11469566063dSJacob Faibussowitsch         PetscCall(VecDuplicate(mglevels[i]->b,&mglevels[i]->crb));
11479566063dSJacob Faibussowitsch         PetscCall(VecZeroEntries(mglevels[i]->crb));
114841b6fd38SMatthew G. Knepley       }
1149f3fbd535SBarry Smith     }
1150f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
1151f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
1152f3fbd535SBarry Smith       Vec *vec;
11539566063dSJacob Faibussowitsch       PetscCall(KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL));
11549566063dSJacob Faibussowitsch       PetscCall(PCMGSetR(pc,n-1,*vec));
11559566063dSJacob Faibussowitsch       PetscCall(VecDestroy(vec));
11569566063dSJacob Faibussowitsch       PetscCall(PetscFree(vec));
1157f3fbd535SBarry Smith     }
115841b6fd38SMatthew G. Knepley     if (doCR) {
11599566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx));
11609566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb));
11619566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(mglevels[n-1]->crb));
116241b6fd38SMatthew G. Knepley     }
1163f3fbd535SBarry Smith   }
1164f3fbd535SBarry Smith 
116598e04984SBarry Smith   if (pc->dm) {
116698e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
116798e04984SBarry Smith     for (i=0; i<n-1; i++) {
116898e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
116998e04984SBarry Smith     }
117098e04984SBarry Smith   }
1171f3fbd535SBarry Smith 
1172f3fbd535SBarry Smith   for (i=1; i<n; i++) {
11732a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1) {
1174f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
11759566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE));
1176f3fbd535SBarry Smith     }
11779566063dSJacob Faibussowitsch     if (mglevels[i]->cr) PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE));
11789566063dSJacob Faibussowitsch     if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0));
11799566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(mglevels[i]->smoothd));
1180c0decd05SBarry Smith     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1181899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1182899639b0SHong Zhang     }
11839566063dSJacob Faibussowitsch     if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0));
1184d42688cbSBarry Smith     if (!mglevels[i]->residual) {
1185d42688cbSBarry Smith       Mat mat;
11869566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothd,&mat,NULL));
11879566063dSJacob Faibussowitsch       PetscCall(PCMGSetResidual(pc,i,PCMGResidualDefault,mat));
1188d42688cbSBarry Smith     }
1189fcb023d4SJed Brown     if (!mglevels[i]->residualtranspose) {
1190fcb023d4SJed Brown       Mat mat;
11919566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(mglevels[i]->smoothd,&mat,NULL));
11929566063dSJacob Faibussowitsch       PetscCall(PCMGSetResidualTranspose(pc,i,PCMGResidualTransposeDefault,mat));
1193fcb023d4SJed Brown     }
1194f3fbd535SBarry Smith   }
1195f3fbd535SBarry Smith   for (i=1; i<n; i++) {
1196f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
1197f3fbd535SBarry Smith       Mat downmat,downpmat;
1198f3fbd535SBarry Smith 
1199f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
12009566063dSJacob Faibussowitsch       PetscCall(KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL));
1201f3fbd535SBarry Smith       if (!opsset) {
12029566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat));
12039566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat));
1204f3fbd535SBarry Smith       }
1205f3fbd535SBarry Smith 
12069566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE));
12079566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0));
12089566063dSJacob Faibussowitsch       PetscCall(KSPSetUp(mglevels[i]->smoothu));
1209c0decd05SBarry Smith       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
1210899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
1211899639b0SHong Zhang       }
12129566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0));
1213f3fbd535SBarry Smith     }
121441b6fd38SMatthew G. Knepley     if (mglevels[i]->cr) {
121541b6fd38SMatthew G. Knepley       Mat downmat,downpmat;
121641b6fd38SMatthew G. Knepley 
121741b6fd38SMatthew G. Knepley       /* check if operators have been set for up, if not use down operators to set them */
12189566063dSJacob Faibussowitsch       PetscCall(KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL));
121941b6fd38SMatthew G. Knepley       if (!opsset) {
12209566063dSJacob Faibussowitsch         PetscCall(KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat));
12219566063dSJacob Faibussowitsch         PetscCall(KSPSetOperators(mglevels[i]->cr,downmat,downpmat));
122241b6fd38SMatthew G. Knepley       }
122341b6fd38SMatthew G. Knepley 
12249566063dSJacob Faibussowitsch       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE));
12259566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0));
12269566063dSJacob Faibussowitsch       PetscCall(KSPSetUp(mglevels[i]->cr));
122741b6fd38SMatthew G. Knepley       if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) {
122841b6fd38SMatthew G. Knepley         pc->failedreason = PC_SUBPC_ERROR;
122941b6fd38SMatthew G. Knepley       }
12309566063dSJacob Faibussowitsch       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0));
123141b6fd38SMatthew G. Knepley     }
1232f3fbd535SBarry Smith   }
1233f3fbd535SBarry Smith 
12349566063dSJacob Faibussowitsch   if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0));
12359566063dSJacob Faibussowitsch   PetscCall(KSPSetUp(mglevels[0]->smoothd));
1236c0decd05SBarry Smith   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1237899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
1238899639b0SHong Zhang   }
12399566063dSJacob Faibussowitsch   if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0));
1240f3fbd535SBarry Smith 
1241f3fbd535SBarry Smith   /*
1242f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
1243e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
1244f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
1245f3fbd535SBarry Smith 
1246f3fbd535SBarry Smith    Only support one or the other at the same time.
1247f3fbd535SBarry Smith   */
1248f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
12499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL));
1250ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
1251f3fbd535SBarry Smith   dump = PETSC_FALSE;
1252f3fbd535SBarry Smith #endif
12539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL));
1254ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
1255f3fbd535SBarry Smith 
1256f3fbd535SBarry Smith   if (viewer) {
1257f3fbd535SBarry Smith     for (i=1; i<n; i++) {
12589566063dSJacob Faibussowitsch       PetscCall(MatView(mglevels[i]->restrct,viewer));
1259f3fbd535SBarry Smith     }
1260f3fbd535SBarry Smith     for (i=0; i<n; i++) {
12619566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(mglevels[i]->smoothd,&pc));
12629566063dSJacob Faibussowitsch       PetscCall(MatView(pc->mat,viewer));
1263f3fbd535SBarry Smith     }
1264f3fbd535SBarry Smith   }
1265f3fbd535SBarry Smith   PetscFunctionReturn(0);
1266f3fbd535SBarry Smith }
1267f3fbd535SBarry Smith 
1268f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
1269f3fbd535SBarry Smith 
127097d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
127197d33e41SMatthew G. Knepley {
127297d33e41SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
127397d33e41SMatthew G. Knepley 
127497d33e41SMatthew G. Knepley   PetscFunctionBegin;
127597d33e41SMatthew G. Knepley   *levels = mg->nlevels;
127697d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
127797d33e41SMatthew G. Knepley }
127897d33e41SMatthew G. Knepley 
12794b9ad928SBarry Smith /*@
128097177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
12814b9ad928SBarry Smith 
12824b9ad928SBarry Smith    Not Collective
12834b9ad928SBarry Smith 
12844b9ad928SBarry Smith    Input Parameter:
12854b9ad928SBarry Smith .  pc - the preconditioner context
12864b9ad928SBarry Smith 
12874b9ad928SBarry Smith    Output parameter:
12884b9ad928SBarry Smith .  levels - the number of levels
12894b9ad928SBarry Smith 
12904b9ad928SBarry Smith    Level: advanced
12914b9ad928SBarry Smith 
129297177400SBarry Smith .seealso: PCMGSetLevels()
12934b9ad928SBarry Smith @*/
12947087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
12954b9ad928SBarry Smith {
12964b9ad928SBarry Smith   PetscFunctionBegin;
12970700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12984482741eSBarry Smith   PetscValidIntPointer(levels,2);
129997d33e41SMatthew G. Knepley   *levels = 0;
13009566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels)));
13014b9ad928SBarry Smith   PetscFunctionReturn(0);
13024b9ad928SBarry Smith }
13034b9ad928SBarry Smith 
13044b9ad928SBarry Smith /*@
1305e7d4b4cbSMark Adams    PCMGGetGridComplexity - compute operator and grid complexity of MG hierarchy
1306e7d4b4cbSMark Adams 
1307e7d4b4cbSMark Adams    Input Parameter:
1308e7d4b4cbSMark Adams .  pc - the preconditioner context
1309e7d4b4cbSMark Adams 
131090db8557SMark Adams    Output Parameters:
131190db8557SMark Adams +  gc - grid complexity = sum_i(n_i) / n_0
131290db8557SMark Adams -  oc - operator complexity = sum_i(nnz_i) / nnz_0
1313e7d4b4cbSMark Adams 
1314e7d4b4cbSMark Adams    Level: advanced
1315e7d4b4cbSMark Adams 
1316e7d4b4cbSMark Adams .seealso: PCMGGetLevels()
1317e7d4b4cbSMark Adams @*/
1318e7d4b4cbSMark Adams PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc, PetscReal *oc)
1319e7d4b4cbSMark Adams {
1320e7d4b4cbSMark Adams   PC_MG          *mg      = (PC_MG*)pc->data;
1321e7d4b4cbSMark Adams   PC_MG_Levels   **mglevels = mg->levels;
1322e7d4b4cbSMark Adams   PetscInt       lev,N;
1323e7d4b4cbSMark Adams   PetscLogDouble nnz0 = 0, sgc = 0, soc = 0, n0 = 0;
1324e7d4b4cbSMark Adams   MatInfo        info;
1325e7d4b4cbSMark Adams 
1326e7d4b4cbSMark Adams   PetscFunctionBegin;
132790db8557SMark Adams   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
132890db8557SMark Adams   if (gc) PetscValidRealPointer(gc,2);
132990db8557SMark Adams   if (oc) PetscValidRealPointer(oc,3);
1330e7d4b4cbSMark Adams   if (!pc->setupcalled) {
133190db8557SMark Adams     if (gc) *gc = 0;
133290db8557SMark Adams     if (oc) *oc = 0;
1333e7d4b4cbSMark Adams     PetscFunctionReturn(0);
1334e7d4b4cbSMark Adams   }
1335e7d4b4cbSMark Adams   PetscCheck(mg->nlevels > 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels");
1336e7d4b4cbSMark Adams   for (lev=0; lev<mg->nlevels; lev++) {
1337e7d4b4cbSMark Adams     Mat dB;
13389566063dSJacob Faibussowitsch     PetscCall(KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB));
13399566063dSJacob Faibussowitsch     PetscCall(MatGetInfo(dB,MAT_GLOBAL_SUM,&info)); /* global reduction */
13409566063dSJacob Faibussowitsch     PetscCall(MatGetSize(dB,&N,NULL));
1341e7d4b4cbSMark Adams     sgc += N;
1342e7d4b4cbSMark Adams     soc += info.nz_used;
1343e7d4b4cbSMark Adams     if (lev==mg->nlevels-1) {nnz0 = info.nz_used; n0 = N;}
1344e7d4b4cbSMark Adams   }
134590db8557SMark Adams   if (n0 > 0 && gc) *gc = (PetscReal)(sgc/n0);
1346e7d4b4cbSMark Adams   else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number for grid points on finest level is not available");
134790db8557SMark Adams   if (nnz0 > 0 && oc) *oc = (PetscReal)(soc/nnz0);
1348e7d4b4cbSMark Adams   PetscFunctionReturn(0);
1349e7d4b4cbSMark Adams }
1350e7d4b4cbSMark Adams 
1351e7d4b4cbSMark Adams /*@
135297177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
13534b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
13544b9ad928SBarry Smith 
1355ad4df100SBarry Smith    Logically Collective on PC
13564b9ad928SBarry Smith 
13574b9ad928SBarry Smith    Input Parameters:
13584b9ad928SBarry Smith +  pc - the preconditioner context
13599dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
13609dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
13614b9ad928SBarry Smith 
13624b9ad928SBarry Smith    Options Database Key:
13634b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
13644b9ad928SBarry Smith    additive, full, kaskade
13654b9ad928SBarry Smith 
13664b9ad928SBarry Smith    Level: advanced
13674b9ad928SBarry Smith 
136897177400SBarry Smith .seealso: PCMGSetLevels()
13694b9ad928SBarry Smith @*/
13707087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
13714b9ad928SBarry Smith {
1372f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
13734b9ad928SBarry Smith 
13744b9ad928SBarry Smith   PetscFunctionBegin;
13750700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1376c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
137731567311SBarry Smith   mg->am = form;
13789dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1379c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
1380c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1381c60c7ad4SBarry Smith }
1382c60c7ad4SBarry Smith 
1383c60c7ad4SBarry Smith /*@
1384c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
1385c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
1386c60c7ad4SBarry Smith 
1387c60c7ad4SBarry Smith    Logically Collective on PC
1388c60c7ad4SBarry Smith 
1389c60c7ad4SBarry Smith    Input Parameter:
1390c60c7ad4SBarry Smith .  pc - the preconditioner context
1391c60c7ad4SBarry Smith 
1392c60c7ad4SBarry Smith    Output Parameter:
1393c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
1394c60c7ad4SBarry Smith 
1395c60c7ad4SBarry Smith    Level: advanced
1396c60c7ad4SBarry Smith 
1397c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
1398c60c7ad4SBarry Smith @*/
1399c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1400c60c7ad4SBarry Smith {
1401c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1402c60c7ad4SBarry Smith 
1403c60c7ad4SBarry Smith   PetscFunctionBegin;
1404c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1405c60c7ad4SBarry Smith   *type = mg->am;
14064b9ad928SBarry Smith   PetscFunctionReturn(0);
14074b9ad928SBarry Smith }
14084b9ad928SBarry Smith 
14094b9ad928SBarry Smith /*@
14100d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
14114b9ad928SBarry Smith    complicated cycling.
14124b9ad928SBarry Smith 
1413ad4df100SBarry Smith    Logically Collective on PC
14144b9ad928SBarry Smith 
14154b9ad928SBarry Smith    Input Parameters:
1416c2be2410SBarry Smith +  pc - the multigrid context
1417c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
14184b9ad928SBarry Smith 
14194b9ad928SBarry Smith    Options Database Key:
1420c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
14214b9ad928SBarry Smith 
14224b9ad928SBarry Smith    Level: advanced
14234b9ad928SBarry Smith 
14240d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
14254b9ad928SBarry Smith @*/
14267087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
14274b9ad928SBarry Smith {
1428f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
1429f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
143079416396SBarry Smith   PetscInt     i,levels;
14314b9ad928SBarry Smith 
14324b9ad928SBarry Smith   PetscFunctionBegin;
14330700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
143469fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
143528b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1436f3fbd535SBarry Smith   levels = mglevels[0]->levels;
14372fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
14384b9ad928SBarry Smith   PetscFunctionReturn(0);
14394b9ad928SBarry Smith }
14404b9ad928SBarry Smith 
14418cc2d5dfSBarry Smith /*@
14428cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
14438cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
14448cc2d5dfSBarry Smith 
1445ad4df100SBarry Smith    Logically Collective on PC
14468cc2d5dfSBarry Smith 
14478cc2d5dfSBarry Smith    Input Parameters:
14488cc2d5dfSBarry Smith +  pc - the multigrid context
14498cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
14508cc2d5dfSBarry Smith 
14518cc2d5dfSBarry Smith    Options Database Key:
1452147403d9SBarry Smith .  -pc_mg_multiplicative_cycles n - set the number of cycles
14538cc2d5dfSBarry Smith 
14548cc2d5dfSBarry Smith    Level: advanced
14558cc2d5dfSBarry Smith 
145695452b02SPatrick Sanan    Notes:
145795452b02SPatrick Sanan     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
14588cc2d5dfSBarry Smith 
14598cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
14608cc2d5dfSBarry Smith @*/
14617087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
14628cc2d5dfSBarry Smith {
1463f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
14648cc2d5dfSBarry Smith 
14658cc2d5dfSBarry Smith   PetscFunctionBegin;
14660700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1467c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
14682a21e185SBarry Smith   mg->cyclesperpcapply = n;
14698cc2d5dfSBarry Smith   PetscFunctionReturn(0);
14708cc2d5dfSBarry Smith }
14718cc2d5dfSBarry Smith 
14722134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1473fb15c04eSBarry Smith {
1474fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1475fb15c04eSBarry Smith 
1476fb15c04eSBarry Smith   PetscFunctionBegin;
14772134b1e4SBarry Smith   mg->galerkin = use;
1478fb15c04eSBarry Smith   PetscFunctionReturn(0);
1479fb15c04eSBarry Smith }
1480fb15c04eSBarry Smith 
1481c2be2410SBarry Smith /*@
148297177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
148382b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1484c2be2410SBarry Smith 
1485ad4df100SBarry Smith    Logically Collective on PC
1486c2be2410SBarry Smith 
1487c2be2410SBarry Smith    Input Parameters:
1488c91913d3SJed Brown +  pc - the multigrid context
14892134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1490c2be2410SBarry Smith 
1491c2be2410SBarry Smith    Options Database Key:
1492147403d9SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - set the matrices to form via the Galerkin process
1493c2be2410SBarry Smith 
1494c2be2410SBarry Smith    Level: intermediate
1495c2be2410SBarry Smith 
149695452b02SPatrick Sanan    Notes:
149795452b02SPatrick Sanan     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
14981c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
14991c1aac46SBarry Smith 
15002134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
15013fc8bf9cSBarry Smith 
1502c2be2410SBarry Smith @*/
15032134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1504c2be2410SBarry Smith {
1505c2be2410SBarry Smith   PetscFunctionBegin;
15060700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
15079566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use)));
1508c2be2410SBarry Smith   PetscFunctionReturn(0);
1509c2be2410SBarry Smith }
1510c2be2410SBarry Smith 
15113fc8bf9cSBarry Smith /*@
15123fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
151382b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
15143fc8bf9cSBarry Smith 
15153fc8bf9cSBarry Smith    Not Collective
15163fc8bf9cSBarry Smith 
15173fc8bf9cSBarry Smith    Input Parameter:
15183fc8bf9cSBarry Smith .  pc - the multigrid context
15193fc8bf9cSBarry Smith 
15203fc8bf9cSBarry Smith    Output Parameter:
15212134b1e4SBarry 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
15223fc8bf9cSBarry Smith 
15233fc8bf9cSBarry Smith    Level: intermediate
15243fc8bf9cSBarry Smith 
15252134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
15263fc8bf9cSBarry Smith 
15273fc8bf9cSBarry Smith @*/
15282134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
15293fc8bf9cSBarry Smith {
1530f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
15313fc8bf9cSBarry Smith 
15323fc8bf9cSBarry Smith   PetscFunctionBegin;
15330700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
15342134b1e4SBarry Smith   *galerkin = mg->galerkin;
15353fc8bf9cSBarry Smith   PetscFunctionReturn(0);
15363fc8bf9cSBarry Smith }
15373fc8bf9cSBarry Smith 
1538f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1539f3b08a26SMatthew G. Knepley {
1540f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
1541f3b08a26SMatthew G. Knepley 
1542f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1543f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = adapt;
1544f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1545f3b08a26SMatthew G. Knepley }
1546f3b08a26SMatthew G. Knepley 
1547f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1548f3b08a26SMatthew G. Knepley {
1549f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
1550f3b08a26SMatthew G. Knepley 
1551f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1552f3b08a26SMatthew G. Knepley   *adapt = mg->adaptInterpolation;
1553f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1554f3b08a26SMatthew G. Knepley }
1555f3b08a26SMatthew G. Knepley 
155641b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr)
155741b6fd38SMatthew G. Knepley {
155841b6fd38SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
155941b6fd38SMatthew G. Knepley 
156041b6fd38SMatthew G. Knepley   PetscFunctionBegin;
156141b6fd38SMatthew G. Knepley   mg->compatibleRelaxation = cr;
156241b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
156341b6fd38SMatthew G. Knepley }
156441b6fd38SMatthew G. Knepley 
156541b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr)
156641b6fd38SMatthew G. Knepley {
156741b6fd38SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
156841b6fd38SMatthew G. Knepley 
156941b6fd38SMatthew G. Knepley   PetscFunctionBegin;
157041b6fd38SMatthew G. Knepley   *cr = mg->compatibleRelaxation;
157141b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
157241b6fd38SMatthew G. Knepley }
157341b6fd38SMatthew G. Knepley 
1574f3b08a26SMatthew G. Knepley /*@
1575f3b08a26SMatthew 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.
1576f3b08a26SMatthew G. Knepley 
1577f3b08a26SMatthew G. Knepley   Logically Collective on PC
1578f3b08a26SMatthew G. Knepley 
1579f3b08a26SMatthew G. Knepley   Input Parameters:
1580f3b08a26SMatthew G. Knepley + pc    - the multigrid context
1581f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator
1582f3b08a26SMatthew G. Knepley 
1583f3b08a26SMatthew G. Knepley   Options Database Keys:
1584f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp                     - Turn on adaptation
1585f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int>             - The number of modes to use, should be divisible by dimension
1586f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector
1587f3b08a26SMatthew G. Knepley 
1588f3b08a26SMatthew G. Knepley   Level: intermediate
1589f3b08a26SMatthew G. Knepley 
1590f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin
1591f3b08a26SMatthew G. Knepley .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1592f3b08a26SMatthew G. Knepley @*/
1593f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1594f3b08a26SMatthew G. Knepley {
1595f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1596f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
15979566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt)));
1598f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1599f3b08a26SMatthew G. Knepley }
1600f3b08a26SMatthew G. Knepley 
1601f3b08a26SMatthew G. Knepley /*@
1602f3b08a26SMatthew G. Knepley   PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated.
1603f3b08a26SMatthew G. Knepley 
1604f3b08a26SMatthew G. Knepley   Not collective
1605f3b08a26SMatthew G. Knepley 
1606f3b08a26SMatthew G. Knepley   Input Parameter:
1607f3b08a26SMatthew G. Knepley . pc    - the multigrid context
1608f3b08a26SMatthew G. Knepley 
1609f3b08a26SMatthew G. Knepley   Output Parameter:
1610f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator
1611f3b08a26SMatthew G. Knepley 
1612f3b08a26SMatthew G. Knepley   Level: intermediate
1613f3b08a26SMatthew G. Knepley 
1614f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin
1615f3b08a26SMatthew G. Knepley .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1616f3b08a26SMatthew G. Knepley @*/
1617f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1618f3b08a26SMatthew G. Knepley {
1619f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1620f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1621dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(adapt, 2);
16229566063dSJacob Faibussowitsch   PetscCall(PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt)));
1623f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1624f3b08a26SMatthew G. Knepley }
1625f3b08a26SMatthew G. Knepley 
16264b9ad928SBarry Smith /*@
162741b6fd38SMatthew G. Knepley   PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation.
162841b6fd38SMatthew G. Knepley 
162941b6fd38SMatthew G. Knepley   Logically Collective on PC
163041b6fd38SMatthew G. Knepley 
163141b6fd38SMatthew G. Knepley   Input Parameters:
163241b6fd38SMatthew G. Knepley + pc - the multigrid context
163341b6fd38SMatthew G. Knepley - cr - flag for compatible relaxation
163441b6fd38SMatthew G. Knepley 
163541b6fd38SMatthew G. Knepley   Options Database Keys:
163641b6fd38SMatthew G. Knepley . -pc_mg_adapt_cr - Turn on compatible relaxation
163741b6fd38SMatthew G. Knepley 
163841b6fd38SMatthew G. Knepley   Level: intermediate
163941b6fd38SMatthew G. Knepley 
164041b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin
164141b6fd38SMatthew G. Knepley .seealso: PCMGGetAdaptCR(), PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
164241b6fd38SMatthew G. Knepley @*/
164341b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr)
164441b6fd38SMatthew G. Knepley {
164541b6fd38SMatthew G. Knepley   PetscFunctionBegin;
164641b6fd38SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16479566063dSJacob Faibussowitsch   PetscCall(PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr)));
164841b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
164941b6fd38SMatthew G. Knepley }
165041b6fd38SMatthew G. Knepley 
165141b6fd38SMatthew G. Knepley /*@
165241b6fd38SMatthew G. Knepley   PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation.
165341b6fd38SMatthew G. Knepley 
165441b6fd38SMatthew G. Knepley   Not collective
165541b6fd38SMatthew G. Knepley 
165641b6fd38SMatthew G. Knepley   Input Parameter:
165741b6fd38SMatthew G. Knepley . pc    - the multigrid context
165841b6fd38SMatthew G. Knepley 
165941b6fd38SMatthew G. Knepley   Output Parameter:
166041b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion
166141b6fd38SMatthew G. Knepley 
166241b6fd38SMatthew G. Knepley   Level: intermediate
166341b6fd38SMatthew G. Knepley 
166441b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin
166541b6fd38SMatthew G. Knepley .seealso: PCMGSetAdaptCR(), PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
166641b6fd38SMatthew G. Knepley @*/
166741b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr)
166841b6fd38SMatthew G. Knepley {
166941b6fd38SMatthew G. Knepley   PetscFunctionBegin;
167041b6fd38SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1671dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(cr, 2);
16729566063dSJacob Faibussowitsch   PetscCall(PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr)));
167341b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
167441b6fd38SMatthew G. Knepley }
167541b6fd38SMatthew G. Knepley 
167641b6fd38SMatthew G. Knepley /*@
167706792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1678710315b6SLawrence Mitchell    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1679710315b6SLawrence Mitchell    pre- and post-smoothing steps.
168006792cafSBarry Smith 
168106792cafSBarry Smith    Logically Collective on PC
168206792cafSBarry Smith 
168306792cafSBarry Smith    Input Parameters:
168406792cafSBarry Smith +  mg - the multigrid context
168506792cafSBarry Smith -  n - the number of smoothing steps
168606792cafSBarry Smith 
168706792cafSBarry Smith    Options Database Key:
1688a2b725a8SWilliam Gropp .  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
168906792cafSBarry Smith 
169006792cafSBarry Smith    Level: advanced
169106792cafSBarry Smith 
169295452b02SPatrick Sanan    Notes:
169395452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
169406792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
169506792cafSBarry Smith 
1696710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp()
169706792cafSBarry Smith @*/
169806792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
169906792cafSBarry Smith {
170006792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
170106792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
170206792cafSBarry Smith   PetscInt       i,levels;
170306792cafSBarry Smith 
170406792cafSBarry Smith   PetscFunctionBegin;
170506792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
170606792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
170728b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
170806792cafSBarry Smith   levels = mglevels[0]->levels;
170906792cafSBarry Smith 
171006792cafSBarry Smith   for (i=1; i<levels; i++) {
17119566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n));
17129566063dSJacob Faibussowitsch     PetscCall(KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n));
171306792cafSBarry Smith     mg->default_smoothu = n;
171406792cafSBarry Smith     mg->default_smoothd = n;
171506792cafSBarry Smith   }
171606792cafSBarry Smith   PetscFunctionReturn(0);
171706792cafSBarry Smith }
171806792cafSBarry Smith 
1719f442ab6aSBarry Smith /*@
17208e5aa403SBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels
1721710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1722f442ab6aSBarry Smith 
1723f442ab6aSBarry Smith    Logically Collective on PC
1724f442ab6aSBarry Smith 
1725f442ab6aSBarry Smith    Input Parameters:
1726f442ab6aSBarry Smith .  pc - the preconditioner context
1727f442ab6aSBarry Smith 
1728f442ab6aSBarry Smith    Options Database Key:
1729147403d9SBarry Smith .  -pc_mg_distinct_smoothup <bool> - use distinct smoothing objects
1730f442ab6aSBarry Smith 
1731f442ab6aSBarry Smith    Level: advanced
1732f442ab6aSBarry Smith 
173395452b02SPatrick Sanan    Notes:
173495452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
1735f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1736f442ab6aSBarry Smith 
1737710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth()
1738f442ab6aSBarry Smith @*/
1739f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1740f442ab6aSBarry Smith {
1741f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1742f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1743f442ab6aSBarry Smith   PetscInt       i,levels;
1744f442ab6aSBarry Smith   KSP            subksp;
1745f442ab6aSBarry Smith 
1746f442ab6aSBarry Smith   PetscFunctionBegin;
1747f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
174828b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1749f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1750f442ab6aSBarry Smith 
1751f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1752710315b6SLawrence Mitchell     const char *prefix = NULL;
1753f442ab6aSBarry Smith     /* make sure smoother up and down are different */
17549566063dSJacob Faibussowitsch     PetscCall(PCMGGetSmootherUp(pc,i,&subksp));
17559566063dSJacob Faibussowitsch     PetscCall(KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix));
17569566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(subksp,prefix));
17579566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(subksp,"up_"));
1758f442ab6aSBarry Smith   }
1759f442ab6aSBarry Smith   PetscFunctionReturn(0);
1760f442ab6aSBarry Smith }
1761f442ab6aSBarry Smith 
176207a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
176307a4832bSFande Kong PetscErrorCode  PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[])
176407a4832bSFande Kong {
176507a4832bSFande Kong   PC_MG          *mg        = (PC_MG*)pc->data;
176607a4832bSFande Kong   PC_MG_Levels   **mglevels = mg->levels;
176707a4832bSFande Kong   Mat            *mat;
176807a4832bSFande Kong   PetscInt       l;
176907a4832bSFande Kong 
177007a4832bSFande Kong   PetscFunctionBegin;
177128b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
17729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mg->nlevels,&mat));
177307a4832bSFande Kong   for (l=1; l< mg->nlevels; l++) {
177407a4832bSFande Kong     mat[l-1] = mglevels[l]->interpolate;
17759566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat[l-1]));
177607a4832bSFande Kong   }
177707a4832bSFande Kong   *num_levels = mg->nlevels;
177807a4832bSFande Kong   *interpolations = mat;
177907a4832bSFande Kong   PetscFunctionReturn(0);
178007a4832bSFande Kong }
178107a4832bSFande Kong 
178207a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
178307a4832bSFande Kong PetscErrorCode  PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[])
178407a4832bSFande Kong {
178507a4832bSFande Kong   PC_MG          *mg        = (PC_MG*)pc->data;
178607a4832bSFande Kong   PC_MG_Levels   **mglevels = mg->levels;
178707a4832bSFande Kong   PetscInt       l;
178807a4832bSFande Kong   Mat            *mat;
178907a4832bSFande Kong 
179007a4832bSFande Kong   PetscFunctionBegin;
179128b400f6SJacob Faibussowitsch   PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
17929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mg->nlevels,&mat));
179307a4832bSFande Kong   for (l=0; l<mg->nlevels-1; l++) {
17949566063dSJacob Faibussowitsch     PetscCall(KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l])));
17959566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)mat[l]));
179607a4832bSFande Kong   }
179707a4832bSFande Kong   *num_levels = mg->nlevels;
179807a4832bSFande Kong   *coarseOperators = mat;
179907a4832bSFande Kong   PetscFunctionReturn(0);
180007a4832bSFande Kong }
180107a4832bSFande Kong 
1802f3b08a26SMatthew G. Knepley /*@C
1803f3b08a26SMatthew G. Knepley   PCMGRegisterCoarseSpaceConstructor -  Adds a method to the PCMG package for coarse space construction.
1804f3b08a26SMatthew G. Knepley 
1805f3b08a26SMatthew G. Knepley   Not collective
1806f3b08a26SMatthew G. Knepley 
1807f3b08a26SMatthew G. Knepley   Input Parameters:
1808f3b08a26SMatthew G. Knepley + name     - name of the constructor
1809f3b08a26SMatthew G. Knepley - function - constructor routine
1810f3b08a26SMatthew G. Knepley 
1811f3b08a26SMatthew G. Knepley   Notes:
1812f3b08a26SMatthew G. Knepley   Calling sequence for the routine:
1813f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1814f3b08a26SMatthew G. Knepley $   pc        - The PC object
1815f3b08a26SMatthew G. Knepley $   l         - The multigrid level, 0 is the coarse level
1816f3b08a26SMatthew G. Knepley $   dm        - The DM for this level
1817f3b08a26SMatthew G. Knepley $   smooth    - The level smoother
1818f3b08a26SMatthew G. Knepley $   Nc        - The size of the coarse space
1819f3b08a26SMatthew G. Knepley $   initGuess - Basis for an initial guess for the space
1820f3b08a26SMatthew G. Knepley $   coarseSp  - A basis for the computed coarse space
1821f3b08a26SMatthew G. Knepley 
1822f3b08a26SMatthew G. Knepley   Level: advanced
1823f3b08a26SMatthew G. Knepley 
1824f3b08a26SMatthew G. Knepley .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister()
1825f3b08a26SMatthew G. Knepley @*/
1826f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1827f3b08a26SMatthew G. Knepley {
1828f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
18299566063dSJacob Faibussowitsch   PetscCall(PCInitializePackage());
18309566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&PCMGCoarseList,name,function));
1831f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1832f3b08a26SMatthew G. Knepley }
1833f3b08a26SMatthew G. Knepley 
1834f3b08a26SMatthew G. Knepley /*@C
1835f3b08a26SMatthew G. Knepley   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.
1836f3b08a26SMatthew G. Knepley 
1837f3b08a26SMatthew G. Knepley   Not collective
1838f3b08a26SMatthew G. Knepley 
1839f3b08a26SMatthew G. Knepley   Input Parameter:
1840f3b08a26SMatthew G. Knepley . name     - name of the constructor
1841f3b08a26SMatthew G. Knepley 
1842f3b08a26SMatthew G. Knepley   Output Parameter:
1843f3b08a26SMatthew G. Knepley . function - constructor routine
1844f3b08a26SMatthew G. Knepley 
1845f3b08a26SMatthew G. Knepley   Notes:
1846f3b08a26SMatthew G. Knepley   Calling sequence for the routine:
1847f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1848f3b08a26SMatthew G. Knepley $   pc        - The PC object
1849f3b08a26SMatthew G. Knepley $   l         - The multigrid level, 0 is the coarse level
1850f3b08a26SMatthew G. Knepley $   dm        - The DM for this level
1851f3b08a26SMatthew G. Knepley $   smooth    - The level smoother
1852f3b08a26SMatthew G. Knepley $   Nc        - The size of the coarse space
1853f3b08a26SMatthew G. Knepley $   initGuess - Basis for an initial guess for the space
1854f3b08a26SMatthew G. Knepley $   coarseSp  - A basis for the computed coarse space
1855f3b08a26SMatthew G. Knepley 
1856f3b08a26SMatthew G. Knepley   Level: advanced
1857f3b08a26SMatthew G. Knepley 
1858f3b08a26SMatthew G. Knepley .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister()
1859f3b08a26SMatthew G. Knepley @*/
1860f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1861f3b08a26SMatthew G. Knepley {
1862f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
18639566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListFind(PCMGCoarseList,name,function));
1864f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1865f3b08a26SMatthew G. Knepley }
1866f3b08a26SMatthew G. Knepley 
18674b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
18684b9ad928SBarry Smith 
18693b09bd56SBarry Smith /*MC
1870ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
18713b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
18723b09bd56SBarry Smith 
18733b09bd56SBarry Smith    Options Database Keys:
18743b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1875391689abSStefano Zampini .  -pc_mg_cycle_type <v,w> - provide the cycle desired
18768c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
18773b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1878710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
18792134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
18808cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
18818cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1882e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
18838cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
18848cf6037eSBarry Smith                         to the binary output file called binaryoutput
18853b09bd56SBarry Smith 
188695452b02SPatrick Sanan    Notes:
18870b3c753eSRichard Tran Mills     If one uses a Krylov method such GMRES or CG as the smoother then one must use KSPFGMRES, KSPGCR, or KSPRICHARDSON as the outer Krylov method
18883b09bd56SBarry Smith 
18898cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
18908cf6037eSBarry Smith 
189123067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
189223067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
189323067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
189423067569SBarry Smith        residual is computed at the end of each cycle.
189523067569SBarry Smith 
18963b09bd56SBarry Smith    Level: intermediate
18973b09bd56SBarry Smith 
18988cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1899710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1900710315b6SLawrence Mitchell            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
190197177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
19020d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
19033b09bd56SBarry Smith M*/
19043b09bd56SBarry Smith 
19058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
19064b9ad928SBarry Smith {
1907f3fbd535SBarry Smith   PC_MG          *mg;
1908f3fbd535SBarry Smith 
19094b9ad928SBarry Smith   PetscFunctionBegin;
19109566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&mg));
19113ec1f749SStefano Zampini   pc->data     = mg;
1912f3fbd535SBarry Smith   mg->nlevels  = -1;
191310eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
19142134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1915f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = PETSC_FALSE;
1916f3b08a26SMatthew G. Knepley   mg->Nc                 = -1;
1917f3b08a26SMatthew G. Knepley   mg->eigenvalue         = -1;
1918f3fbd535SBarry Smith 
191937a44384SMark Adams   pc->useAmat = PETSC_TRUE;
192037a44384SMark Adams 
19214b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
1922fcb023d4SJed Brown   pc->ops->applytranspose = PCApplyTranspose_MG;
192330b0564aSStefano Zampini   pc->ops->matapply       = PCMatApply_MG;
19244b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1925a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
19264b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
19274b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
19284b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1929fb15c04eSBarry Smith 
19309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataRegister(&mg->eigenvalue));
19319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG));
19329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG));
19339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG));
19349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG));
19359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG));
19369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG));
19379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG));
19389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG));
19399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG));
19404b9ad928SBarry Smith   PetscFunctionReturn(0);
19414b9ad928SBarry Smith }
1942