xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision ffc4695bcb29f4b022f59a5fd6bc99fc280ff6d8)
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 
1531567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
164b9ad928SBarry Smith {
1731567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1831567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
196849ba73SBarry Smith   PetscErrorCode ierr;
2031567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
214b9ad928SBarry Smith 
224b9ad928SBarry Smith   PetscFunctionBegin;
2363e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2431567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
25c0decd05SBarry Smith   ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr);
2663e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2731567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2863e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2931567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
3063e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
314b9ad928SBarry Smith 
324b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
3331567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
344b9ad928SBarry Smith       PetscReal rnorm;
3531567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
364b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3770441072SBarry Smith         if (rnorm < mg->abstol) {
384d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
3957622a8eSBarry Smith           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr);
404b9ad928SBarry Smith         } else {
414d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
4257622a8eSBarry Smith           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);CHKERRQ(ierr);
434b9ad928SBarry Smith         }
444b9ad928SBarry Smith         PetscFunctionReturn(0);
454b9ad928SBarry Smith       }
464b9ad928SBarry Smith     }
474b9ad928SBarry Smith 
4831567311SBarry Smith     mgc = *(mglevelsin - 1);
4963e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5031567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
5163e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
52efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
534b9ad928SBarry Smith     while (cycles--) {
5431567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
554b9ad928SBarry Smith     }
5663e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5731567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5863e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5963e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
6031567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
61c0decd05SBarry Smith     ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr);
6241b6fd38SMatthew G. Knepley     if (mglevels->cr) {
6341b6fd38SMatthew G. Knepley       /* TODO Turn on copy and turn off noisy if we have an exact solution
6441b6fd38SMatthew G. Knepley       ierr = VecCopy(mglevels->x, mglevels->crx);CHKERRQ(ierr);
6541b6fd38SMatthew G. Knepley       ierr = VecCopy(mglevels->b, mglevels->crb);CHKERRQ(ierr); */
6641b6fd38SMatthew G. Knepley       ierr = KSPSetNoisy_Private(mglevels->crx);CHKERRQ(ierr);
6741b6fd38SMatthew G. Knepley       ierr = KSPSolve(mglevels->cr,mglevels->crb,mglevels->crx);CHKERRQ(ierr);    /* compatible relaxation */
6841b6fd38SMatthew G. Knepley       ierr = KSPCheckSolve(mglevels->cr,pc,mglevels->crx);CHKERRQ(ierr);
6941b6fd38SMatthew G. Knepley     }
7063e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
714b9ad928SBarry Smith   }
724b9ad928SBarry Smith   PetscFunctionReturn(0);
734b9ad928SBarry Smith }
744b9ad928SBarry Smith 
75ace3abfcSBarry 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)
764b9ad928SBarry Smith {
77f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
78f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
79dfbe8321SBarry Smith   PetscErrorCode ierr;
80391689abSStefano Zampini   PC             tpc;
81391689abSStefano Zampini   PetscBool      changeu,changed;
82f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
834b9ad928SBarry Smith 
844b9ad928SBarry Smith   PetscFunctionBegin;
85a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
86a762d673SBarry Smith   for (i=0; i<levels; i++) {
87a762d673SBarry Smith     if (!mglevels[i]->A) {
88a762d673SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
89a762d673SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
90a762d673SBarry Smith     }
91a762d673SBarry Smith   }
92391689abSStefano Zampini 
93391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
94391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
95391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
96391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
97391689abSStefano Zampini   if (!changed && !changeu) {
98391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
99f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
100391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
101391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
102391689abSStefano Zampini       Vec *vec;
103391689abSStefano Zampini 
104391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
105391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
1067ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
107391689abSStefano Zampini     }
108391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
109391689abSStefano Zampini   }
110f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
1114b9ad928SBarry Smith 
11231567311SBarry Smith   mg->rtol   = rtol;
11331567311SBarry Smith   mg->abstol = abstol;
11431567311SBarry Smith   mg->dtol   = dtol;
1154b9ad928SBarry Smith   if (rtol) {
1164b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
1174b9ad928SBarry Smith     PetscReal rnorm;
1187319c654SBarry Smith     if (zeroguess) {
1197319c654SBarry Smith       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
1207319c654SBarry Smith     } else {
121f3fbd535SBarry Smith       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
1224b9ad928SBarry Smith       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
1237319c654SBarry Smith     }
12431567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
1252fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
1262fa5cd67SKarl Rupp   else mg->ttol = 0.0;
1274b9ad928SBarry Smith 
1284d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
129336babb1SJed Brown      stop prematurely due to small residual */
1304d0a8057SBarry Smith   for (i=1; i<levels; i++) {
131f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
132f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
13323067569SBarry Smith       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
13423067569SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
135f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1364b9ad928SBarry Smith     }
1374d0a8057SBarry Smith   }
1384d0a8057SBarry Smith 
1394d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1404d0a8057SBarry Smith   for (i=0; i<its; i++) {
141f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1424d0a8057SBarry Smith     if (*reason) break;
1434d0a8057SBarry Smith   }
1444d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1454d0a8057SBarry Smith   *outits = i;
146391689abSStefano Zampini   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
1474b9ad928SBarry Smith   PetscFunctionReturn(0);
1484b9ad928SBarry Smith }
1494b9ad928SBarry Smith 
1503aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1513aeaf226SBarry Smith {
1523aeaf226SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1533aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1543aeaf226SBarry Smith   PetscErrorCode ierr;
155f3b08a26SMatthew G. Knepley   PetscInt       i,c,n;
1563aeaf226SBarry Smith 
1573aeaf226SBarry Smith   PetscFunctionBegin;
1583aeaf226SBarry Smith   if (mglevels) {
1593aeaf226SBarry Smith     n = mglevels[0]->levels;
1603aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1613aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1623aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1633aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
16441b6fd38SMatthew G. Knepley       ierr = VecDestroy(&mglevels[i]->crx);CHKERRQ(ierr);
16541b6fd38SMatthew G. Knepley       ierr = VecDestroy(&mglevels[i]->crb);CHKERRQ(ierr);
1663aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1673aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
1680de8493bSLawrence Mitchell       ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr);
16973250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1703aeaf226SBarry Smith     }
17141b6fd38SMatthew G. Knepley     ierr = VecDestroy(&mglevels[n-1]->crx);CHKERRQ(ierr);
17241b6fd38SMatthew G. Knepley     ierr = VecDestroy(&mglevels[n-1]->crb);CHKERRQ(ierr);
173391689abSStefano Zampini     /* this is not null only if the smoother on the finest level
174391689abSStefano Zampini        changes the rhs during PreSolve */
175391689abSStefano Zampini     ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr);
1763aeaf226SBarry Smith 
1773aeaf226SBarry Smith     for (i=0; i<n; i++) {
178f3b08a26SMatthew G. Knepley       if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {ierr = VecDestroy(&mglevels[i]->coarseSpace[c]);CHKERRQ(ierr);}
179f3b08a26SMatthew G. Knepley       ierr = PetscFree(mglevels[i]->coarseSpace);CHKERRQ(ierr);
180f3b08a26SMatthew G. Knepley       mglevels[i]->coarseSpace = NULL;
1813aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1823aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1833aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1843aeaf226SBarry Smith       }
1853aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
18641b6fd38SMatthew G. Knepley       if (mglevels[i]->cr) {ierr = KSPReset(mglevels[i]->cr);CHKERRQ(ierr);}
1873aeaf226SBarry Smith     }
188f3b08a26SMatthew G. Knepley     mg->Nc = 0;
1893aeaf226SBarry Smith   }
1903aeaf226SBarry Smith   PetscFunctionReturn(0);
1913aeaf226SBarry Smith }
1923aeaf226SBarry Smith 
19341b6fd38SMatthew G. Knepley /* Implementing CR
19441b6fd38SMatthew G. Knepley 
19541b6fd38SMatthew 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
19641b6fd38SMatthew G. Knepley 
19741b6fd38SMatthew G. Knepley   Inj^T (Inj Inj^T)^{-1} Inj
19841b6fd38SMatthew G. Knepley 
19941b6fd38SMatthew G. Knepley and if Inj is a VecScatter, as it is now in PETSc, we have
20041b6fd38SMatthew G. Knepley 
20141b6fd38SMatthew G. Knepley   Inj^T Inj
20241b6fd38SMatthew G. Knepley 
20341b6fd38SMatthew G. Knepley and
20441b6fd38SMatthew G. Knepley 
20541b6fd38SMatthew G. Knepley   S = I - Inj^T Inj
20641b6fd38SMatthew G. Knepley 
20741b6fd38SMatthew G. Knepley since
20841b6fd38SMatthew G. Knepley 
20941b6fd38SMatthew G. Knepley   Inj S = Inj - (Inj Inj^T) Inj = 0.
21041b6fd38SMatthew G. Knepley 
21141b6fd38SMatthew G. Knepley Brannick suggests
21241b6fd38SMatthew G. Knepley 
21341b6fd38SMatthew G. Knepley   A \to S^T A S  \qquad\mathrm{and}\qquad M \to S^T M S
21441b6fd38SMatthew G. Knepley 
21541b6fd38SMatthew 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
21641b6fd38SMatthew G. Knepley 
21741b6fd38SMatthew G. Knepley   M^{-1} A \to S M^{-1} A S
21841b6fd38SMatthew G. Knepley 
21941b6fd38SMatthew G. Knepley In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left.
22041b6fd38SMatthew G. Knepley 
22141b6fd38SMatthew G. Knepley   Check: || Inj P - I ||_F < tol
22241b6fd38SMatthew G. Knepley   Check: In general, Inj Inj^T = I
22341b6fd38SMatthew G. Knepley */
22441b6fd38SMatthew G. Knepley 
22541b6fd38SMatthew G. Knepley typedef struct {
22641b6fd38SMatthew G. Knepley   PC       mg;  /* The PCMG object */
22741b6fd38SMatthew G. Knepley   PetscInt l;   /* The multigrid level for this solver */
22841b6fd38SMatthew G. Knepley   Mat      Inj; /* The injection matrix */
22941b6fd38SMatthew G. Knepley   Mat      S;   /* I - Inj^T Inj */
23041b6fd38SMatthew G. Knepley } CRContext;
23141b6fd38SMatthew G. Knepley 
23241b6fd38SMatthew G. Knepley static PetscErrorCode CRSetup_Private(PC pc)
23341b6fd38SMatthew G. Knepley {
23441b6fd38SMatthew G. Knepley   CRContext     *ctx;
23541b6fd38SMatthew G. Knepley   Mat            It;
23641b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
23741b6fd38SMatthew G. Knepley 
23841b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
23941b6fd38SMatthew G. Knepley   ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr);
24041b6fd38SMatthew G. Knepley   ierr = PCMGGetInjection(ctx->mg, ctx->l, &It);CHKERRQ(ierr);
24141b6fd38SMatthew G. Knepley   if (!It) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG");
24241b6fd38SMatthew G. Knepley   ierr = MatCreateTranspose(It, &ctx->Inj);CHKERRQ(ierr);
24341b6fd38SMatthew G. Knepley   ierr = MatCreateNormal(ctx->Inj, &ctx->S);CHKERRQ(ierr);
24441b6fd38SMatthew G. Knepley   ierr = MatScale(ctx->S, -1.0);CHKERRQ(ierr);
24541b6fd38SMatthew G. Knepley   ierr = MatShift(ctx->S,  1.0);CHKERRQ(ierr);
24641b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
24741b6fd38SMatthew G. Knepley }
24841b6fd38SMatthew G. Knepley 
24941b6fd38SMatthew G. Knepley static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y)
25041b6fd38SMatthew G. Knepley {
25141b6fd38SMatthew G. Knepley   CRContext     *ctx;
25241b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
25341b6fd38SMatthew G. Knepley 
25441b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
25541b6fd38SMatthew G. Knepley   ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr);
25641b6fd38SMatthew G. Knepley   ierr = MatMult(ctx->S, x, y);CHKERRQ(ierr);
25741b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
25841b6fd38SMatthew G. Knepley }
25941b6fd38SMatthew G. Knepley 
26041b6fd38SMatthew G. Knepley static PetscErrorCode CRDestroy_Private(PC pc)
26141b6fd38SMatthew G. Knepley {
26241b6fd38SMatthew G. Knepley   CRContext     *ctx;
26341b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
26441b6fd38SMatthew G. Knepley 
26541b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
26641b6fd38SMatthew G. Knepley   ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr);
26741b6fd38SMatthew G. Knepley   ierr = MatDestroy(&ctx->Inj);CHKERRQ(ierr);
26841b6fd38SMatthew G. Knepley   ierr = MatDestroy(&ctx->S);CHKERRQ(ierr);
26941b6fd38SMatthew G. Knepley   ierr = PetscFree(ctx);CHKERRQ(ierr);
27041b6fd38SMatthew G. Knepley   ierr = PCShellSetContext(pc, NULL);CHKERRQ(ierr);
27141b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
27241b6fd38SMatthew G. Knepley }
27341b6fd38SMatthew G. Knepley 
27441b6fd38SMatthew G. Knepley static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr)
27541b6fd38SMatthew G. Knepley {
27641b6fd38SMatthew G. Knepley   CRContext     *ctx;
27741b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
27841b6fd38SMatthew G. Knepley 
27941b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
28041b6fd38SMatthew G. Knepley   ierr = PCCreate(PetscObjectComm((PetscObject) pc), cr);CHKERRQ(ierr);
28141b6fd38SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *cr, "S (complementary projector to injection)");CHKERRQ(ierr);
28241b6fd38SMatthew G. Knepley   ierr = PetscCalloc1(1, &ctx);CHKERRQ(ierr);
28341b6fd38SMatthew G. Knepley   ctx->mg = pc;
28441b6fd38SMatthew G. Knepley   ctx->l  = l;
28541b6fd38SMatthew G. Knepley   ierr = PCSetType(*cr, PCSHELL);CHKERRQ(ierr);
28641b6fd38SMatthew G. Knepley   ierr = PCShellSetContext(*cr, ctx);CHKERRQ(ierr);
28741b6fd38SMatthew G. Knepley   ierr = PCShellSetApply(*cr, CRApply_Private);CHKERRQ(ierr);
28841b6fd38SMatthew G. Knepley   ierr = PCShellSetSetUp(*cr, CRSetup_Private);CHKERRQ(ierr);
28941b6fd38SMatthew G. Knepley   ierr = PCShellSetDestroy(*cr, CRDestroy_Private);CHKERRQ(ierr);
29041b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
29141b6fd38SMatthew G. Knepley }
29241b6fd38SMatthew G. Knepley 
29397d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
2944b9ad928SBarry Smith {
295dfbe8321SBarry Smith   PetscErrorCode ierr;
296f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
297ce94432eSBarry Smith   MPI_Comm       comm;
2983aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
29910eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
30010167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
301f3fbd535SBarry Smith   PetscInt       i;
302f3fbd535SBarry Smith   PetscMPIInt    size;
303f3fbd535SBarry Smith   const char     *prefix;
304f3fbd535SBarry Smith   PC             ipc;
3053aeaf226SBarry Smith   PetscInt       n;
3064b9ad928SBarry Smith 
3074b9ad928SBarry Smith   PetscFunctionBegin;
3080700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
309548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
310548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
3111a97d7b6SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
3123aeaf226SBarry Smith   if (mglevels) {
31310eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
3143aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
3153aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
3163aeaf226SBarry Smith     n    = mglevels[0]->levels;
3173aeaf226SBarry Smith     for (i=0; i<n; i++) {
3183aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
3193aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
3203aeaf226SBarry Smith       }
3213aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
32241b6fd38SMatthew G. Knepley       ierr = KSPDestroy(&mglevels[i]->cr);CHKERRQ(ierr);
3233aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
3243aeaf226SBarry Smith     }
3253aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
3263aeaf226SBarry Smith   }
327f3fbd535SBarry Smith 
328f3fbd535SBarry Smith   mg->nlevels = levels;
329f3fbd535SBarry Smith 
330785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
3313bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
332f3fbd535SBarry Smith 
333f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
334f3fbd535SBarry Smith 
335a9db87a2SMatthew G Knepley   mg->stageApply = 0;
336f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
337b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
3382fa5cd67SKarl Rupp 
339f3fbd535SBarry Smith     mglevels[i]->level               = i;
340f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
34110eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
342336babb1SJed Brown     mg->default_smoothu              = 2;
343336babb1SJed Brown     mg->default_smoothd              = 2;
34463e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
34563e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
34663e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
34763e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
348f3fbd535SBarry Smith 
349f3fbd535SBarry Smith     if (comms) comm = comms[i];
350d5a8f7e6SBarry Smith     if (comm != MPI_COMM_NULL) {
351f3fbd535SBarry Smith       ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
352422a814eSBarry Smith       ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
3530ee9a668SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
3540ee9a668SBarry Smith       ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
3550ee9a668SBarry Smith       ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
3560ee9a668SBarry Smith       if (i || levels == 1) {
3570ee9a668SBarry Smith         char tprefix[128];
3580ee9a668SBarry Smith 
359336babb1SJed Brown         ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
3600059c7bdSJed Brown         ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
361336babb1SJed Brown         ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
362336babb1SJed Brown         ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
363336babb1SJed Brown         ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
3640ee9a668SBarry Smith         ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
365f3fbd535SBarry Smith 
36641b6fd38SMatthew G. Knepley         ierr = PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i);CHKERRQ(ierr);
3670ee9a668SBarry Smith         ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
3680ee9a668SBarry Smith       } else {
369f3fbd535SBarry Smith         ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
370f3fbd535SBarry Smith 
3719dbfc187SHong Zhang         /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
372f3fbd535SBarry Smith         ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
373f3fbd535SBarry Smith         ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
374*ffc4695bSBarry Smith         ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
375f3fbd535SBarry Smith         if (size > 1) {
376f3fbd535SBarry Smith           ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
377f3fbd535SBarry Smith         } else {
378f3fbd535SBarry Smith           ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
379f3fbd535SBarry Smith         }
380753b7fb9SBarry Smith         ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
381f3fbd535SBarry Smith       }
3823bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
383d5a8f7e6SBarry Smith     }
384f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
38531567311SBarry Smith     mg->rtol             = 0.0;
38631567311SBarry Smith     mg->abstol           = 0.0;
38731567311SBarry Smith     mg->dtol             = 0.0;
38831567311SBarry Smith     mg->ttol             = 0.0;
38931567311SBarry Smith     mg->cyclesperpcapply = 1;
390f3fbd535SBarry Smith   }
391f3fbd535SBarry Smith   mg->levels = mglevels;
39210eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
3934b9ad928SBarry Smith   PetscFunctionReturn(0);
3944b9ad928SBarry Smith }
3954b9ad928SBarry Smith 
39697d33e41SMatthew G. Knepley /*@C
39797d33e41SMatthew G. Knepley    PCMGSetLevels - Sets the number of levels to use with MG.
39897d33e41SMatthew G. Knepley    Must be called before any other MG routine.
39997d33e41SMatthew G. Knepley 
40097d33e41SMatthew G. Knepley    Logically Collective on PC
40197d33e41SMatthew G. Knepley 
40297d33e41SMatthew G. Knepley    Input Parameters:
40397d33e41SMatthew G. Knepley +  pc - the preconditioner context
40497d33e41SMatthew G. Knepley .  levels - the number of levels
40597d33e41SMatthew G. Knepley -  comms - optional communicators for each level; this is to allow solving the coarser problems
406d5a8f7e6SBarry Smith            on smaller sets of processes. For processes that are not included in the computation
407d5a8f7e6SBarry Smith            you must pass MPI_COMM_NULL.
40897d33e41SMatthew G. Knepley 
40997d33e41SMatthew G. Knepley    Level: intermediate
41097d33e41SMatthew G. Knepley 
41197d33e41SMatthew G. Knepley    Notes:
41297d33e41SMatthew G. Knepley      If the number of levels is one then the multigrid uses the -mg_levels prefix
41397d33e41SMatthew G. Knepley      for setting the level options rather than the -mg_coarse prefix.
41497d33e41SMatthew G. Knepley 
415d5a8f7e6SBarry Smith      You can free the information in comms after this routine is called.
416d5a8f7e6SBarry Smith 
417d5a8f7e6SBarry Smith      The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level
418d5a8f7e6SBarry Smith      are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on
419d5a8f7e6SBarry Smith      the two levels, rank 0 in the original communicator will pass in an array of 2 communicators
420d5a8f7e6SBarry Smith      of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators
421d5a8f7e6SBarry Smith      the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate
422d5a8f7e6SBarry Smith      in the coarse grid solve.
423d5a8f7e6SBarry Smith 
424d5a8f7e6SBarry Smith      Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one
425d5a8f7e6SBarry Smith      must take special care in providing the restriction and interpolation operation. We recommend
426d5a8f7e6SBarry Smith      providing these as two step operations; first perform a standard restriction or interpolation on
427d5a8f7e6SBarry Smith      the full number of ranks for that level and then use an MPI call to copy the resulting vector
428d5a8f7e6SBarry Smith      array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, not in both
429d5a8f7e6SBarry Smith      cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and
430d5a8f7e6SBarry Smith      recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries.
431d5a8f7e6SBarry Smith 
432d5a8f7e6SBarry Smith 
43397d33e41SMatthew G. Knepley .seealso: PCMGSetType(), PCMGGetLevels()
43497d33e41SMatthew G. Knepley @*/
43597d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
43697d33e41SMatthew G. Knepley {
43797d33e41SMatthew G. Knepley   PetscErrorCode ierr;
43897d33e41SMatthew G. Knepley 
43997d33e41SMatthew G. Knepley   PetscFunctionBegin;
44097d33e41SMatthew G. Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
44197d33e41SMatthew G. Knepley   if (comms) PetscValidPointer(comms,3);
44237b5128cSMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr);
44397d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
44497d33e41SMatthew G. Knepley }
44597d33e41SMatthew G. Knepley 
446c07bf074SBarry Smith 
447c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
448c07bf074SBarry Smith {
449c07bf074SBarry Smith   PetscErrorCode ierr;
450a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
451a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
452a06653b4SBarry Smith   PetscInt       i,n;
453c07bf074SBarry Smith 
454c07bf074SBarry Smith   PetscFunctionBegin;
455a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
456a06653b4SBarry Smith   if (mglevels) {
457a06653b4SBarry Smith     n = mglevels[0]->levels;
458a06653b4SBarry Smith     for (i=0; i<n; i++) {
459a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
4606bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
461a06653b4SBarry Smith       }
4626bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
46341b6fd38SMatthew G. Knepley       ierr = KSPDestroy(&mglevels[i]->cr);CHKERRQ(ierr);
464a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
465a06653b4SBarry Smith     }
466a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
467a06653b4SBarry Smith   }
468c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
469fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);CHKERRQ(ierr);
470fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);CHKERRQ(ierr);
471f3fbd535SBarry Smith   PetscFunctionReturn(0);
472f3fbd535SBarry Smith }
473f3fbd535SBarry Smith 
474f3fbd535SBarry Smith 
475f3fbd535SBarry Smith 
47609573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
47709573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
47809573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
479f3fbd535SBarry Smith 
480f3fbd535SBarry Smith /*
481f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
482f3fbd535SBarry Smith              or full cycle of multigrid.
483f3fbd535SBarry Smith 
484f3fbd535SBarry Smith   Note:
485f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
486f3fbd535SBarry Smith */
487f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
488f3fbd535SBarry Smith {
489f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
490f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
491f3fbd535SBarry Smith   PetscErrorCode ierr;
492391689abSStefano Zampini   PC             tpc;
493f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
494391689abSStefano Zampini   PetscBool      changeu,changed;
495f3fbd535SBarry Smith 
496f3fbd535SBarry Smith   PetscFunctionBegin;
497a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
498e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
499298cc208SBarry Smith   for (i=0; i<levels; i++) {
500e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
50123ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
502298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
503e1d8e5deSBarry Smith     }
504e1d8e5deSBarry Smith   }
505e1d8e5deSBarry Smith 
506391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
507391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
508391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
509391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
510391689abSStefano Zampini   if (!changeu && !changed) {
511391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
512f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
513391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
514391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
515391689abSStefano Zampini       Vec *vec;
516391689abSStefano Zampini 
517391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
518391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
5197ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
520391689abSStefano Zampini     }
521391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
522391689abSStefano Zampini   }
523f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
524391689abSStefano Zampini 
52531567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
526f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
52731567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
5280298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
529f3fbd535SBarry Smith     }
5302fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
53131567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
5322fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
53331567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
5342fa5cd67SKarl Rupp   } else {
535f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
536f3fbd535SBarry Smith   }
537a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
538391689abSStefano Zampini   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
539f3fbd535SBarry Smith   PetscFunctionReturn(0);
540f3fbd535SBarry Smith }
541f3fbd535SBarry Smith 
542f3fbd535SBarry Smith 
5434416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
544f3fbd535SBarry Smith {
545f3fbd535SBarry Smith   PetscErrorCode   ierr;
546710315b6SLawrence Mitchell   PetscInt         levels,cycles;
547f3b08a26SMatthew G. Knepley   PetscBool        flg, flg2;
548f3fbd535SBarry Smith   PC_MG            *mg = (PC_MG*)pc->data;
5493d3eaba7SBarry Smith   PC_MG_Levels     **mglevels;
550f3fbd535SBarry Smith   PCMGType         mgtype;
551f3fbd535SBarry Smith   PCMGCycleType    mgctype;
5522134b1e4SBarry Smith   PCMGGalerkinType gtype;
553f3fbd535SBarry Smith 
554f3fbd535SBarry Smith   PetscFunctionBegin;
5551d743356SStefano Zampini   levels = PetscMax(mg->nlevels,1);
556e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
557f3fbd535SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
5581a97d7b6SStefano Zampini   if (!flg && !mg->levels && pc->dm) {
559eb3f98d2SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
560eb3f98d2SBarry Smith     levels++;
5613aeaf226SBarry Smith     mg->usedmfornumberoflevels = PETSC_TRUE;
562eb3f98d2SBarry Smith   }
5630298fd71SBarry Smith   ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5643aeaf226SBarry Smith   mglevels = mg->levels;
5653aeaf226SBarry Smith 
566f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
567f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
568f3fbd535SBarry Smith   if (flg) {
569f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
5702fa5cd67SKarl Rupp   }
5712134b1e4SBarry Smith   gtype = mg->galerkin;
5722134b1e4SBarry Smith   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
5732134b1e4SBarry Smith   if (flg) {
5742134b1e4SBarry Smith     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
575f3fbd535SBarry Smith   }
576f3b08a26SMatthew G. Knepley   flg2 = PETSC_FALSE;
577f3b08a26SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr);
578f3b08a26SMatthew G. Knepley   if (flg) {ierr = PCMGSetAdaptInterpolation(pc, flg2);CHKERRQ(ierr);}
579f3b08a26SMatthew G. Knepley   ierr = PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg);CHKERRQ(ierr);
580f3b08a26SMatthew G. Knepley   ierr = 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);CHKERRQ(ierr);
581f3b08a26SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);CHKERRQ(ierr);
58241b6fd38SMatthew G. Knepley   flg2 = PETSC_FALSE;
58341b6fd38SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr);
58441b6fd38SMatthew G. Knepley   if (flg) {ierr = PCMGSetAdaptCR(pc, flg2);CHKERRQ(ierr);}
585f56b1265SBarry Smith   flg = PETSC_FALSE;
5868e5aa403SBarry Smith   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
587f442ab6aSBarry Smith   if (flg) {
588f442ab6aSBarry Smith     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
589f442ab6aSBarry Smith   }
59031567311SBarry Smith   mgtype = mg->am;
591f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
592f3fbd535SBarry Smith   if (flg) {
593f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
594f3fbd535SBarry Smith   }
59531567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
5965363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
597f3fbd535SBarry Smith     if (flg) {
598f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
599f3fbd535SBarry Smith     }
600f3fbd535SBarry Smith   }
601f3fbd535SBarry Smith   flg  = PETSC_FALSE;
6020298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
603f3fbd535SBarry Smith   if (flg) {
604f3fbd535SBarry Smith     PetscInt i;
605f3fbd535SBarry Smith     char     eventname[128];
6061a97d7b6SStefano Zampini 
607f3fbd535SBarry Smith     levels = mglevels[0]->levels;
608f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
609f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
61063e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
611f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
61263e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
613f3fbd535SBarry Smith       if (i) {
614f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
61563e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
616f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
61763e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
618f3fbd535SBarry Smith       }
619f3fbd535SBarry Smith     }
620eec35531SMatthew G Knepley 
621ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
622eec35531SMatthew G Knepley     {
623eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
624eec35531SMatthew G Knepley       PetscStageLog stageLog;
625eec35531SMatthew G Knepley       PetscInt      st;
626eec35531SMatthew G Knepley 
627eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
628eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
629eec35531SMatthew G Knepley         PetscBool same;
630eec35531SMatthew G Knepley 
631eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
6322fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
633eec35531SMatthew G Knepley       }
634eec35531SMatthew G Knepley       if (!mg->stageApply) {
635eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
636eec35531SMatthew G Knepley       }
637eec35531SMatthew G Knepley     }
638ec60ca73SSatish Balay #endif
639f3fbd535SBarry Smith   }
640f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
641f3b08a26SMatthew G. Knepley   /* Check option consistency */
642f3b08a26SMatthew G. Knepley   ierr = PCMGGetGalerkin(pc, &gtype);CHKERRQ(ierr);
643f3b08a26SMatthew G. Knepley   ierr = PCMGGetAdaptInterpolation(pc, &flg);CHKERRQ(ierr);
644f3b08a26SMatthew G. Knepley   if (flg && (gtype >= PC_MG_GALERKIN_NONE)) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
645f3fbd535SBarry Smith   PetscFunctionReturn(0);
646f3fbd535SBarry Smith }
647f3fbd535SBarry Smith 
6480a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL};
6490a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL};
6500a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL};
651f3b08a26SMatthew G. Knepley const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL};
652f3fbd535SBarry Smith 
6539804daf3SBarry Smith #include <petscdraw.h>
654f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
655f3fbd535SBarry Smith {
656f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
657f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
658f3fbd535SBarry Smith   PetscErrorCode ierr;
659e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
660219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
661f3fbd535SBarry Smith 
662f3fbd535SBarry Smith   PetscFunctionBegin;
663251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6645b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
665219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
666f3fbd535SBarry Smith   if (iascii) {
667e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
668efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
66931567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
67031567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
671f3fbd535SBarry Smith     }
6722134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
673f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
6742134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
6752134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
6762134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
6772134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
6782134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
6792134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
6804f66f45eSBarry Smith     } else {
6814f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
682f3fbd535SBarry Smith     }
6835adeb434SBarry Smith     if (mg->view){
6845adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
6855adeb434SBarry Smith     }
686f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
687f3fbd535SBarry Smith       if (!i) {
68889cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
689f3fbd535SBarry Smith       } else {
69089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
691f3fbd535SBarry Smith       }
692f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
693f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
694f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
695f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
696f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
697f3fbd535SBarry Smith       } else if (i) {
69889cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
699f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
700f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
701f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
702f3fbd535SBarry Smith       }
70341b6fd38SMatthew G. Knepley       if (i && mglevels[i]->cr) {
70441b6fd38SMatthew G. Knepley         ierr = PetscViewerASCIIPrintf(viewer,"CR solver on level %D -------------------------------\n",i);CHKERRQ(ierr);
70541b6fd38SMatthew G. Knepley         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
70641b6fd38SMatthew G. Knepley         ierr = KSPView(mglevels[i]->cr,viewer);CHKERRQ(ierr);
70741b6fd38SMatthew G. Knepley         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
70841b6fd38SMatthew G. Knepley       }
709f3fbd535SBarry Smith     }
7105b0b0462SBarry Smith   } else if (isbinary) {
7115b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
7125b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
7135b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
7145b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
7155b0b0462SBarry Smith       }
7165b0b0462SBarry Smith     }
717219b1045SBarry Smith   } else if (isdraw) {
718219b1045SBarry Smith     PetscDraw draw;
719b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
720219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
721219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
7220298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
723219b1045SBarry Smith     bottom = y - th;
724219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
725b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
726219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
727219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
728219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
729b4375e8dSPeter Brune       } else {
730b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
731b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
732b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
733b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
734b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
735b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
736b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
737b4375e8dSPeter Brune       }
7380298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
7391cd381d6SBarry Smith       bottom -= th;
740219b1045SBarry Smith     }
7415b0b0462SBarry Smith   }
742f3fbd535SBarry Smith   PetscFunctionReturn(0);
743f3fbd535SBarry Smith }
744f3fbd535SBarry Smith 
745af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
746cab2e9ccSBarry Smith 
747f3fbd535SBarry Smith /*
748f3fbd535SBarry Smith     Calls setup for the KSP on each level
749f3fbd535SBarry Smith */
750f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
751f3fbd535SBarry Smith {
752f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
753f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
754f3fbd535SBarry Smith   PetscErrorCode ierr;
7551a97d7b6SStefano Zampini   PetscInt       i,n;
75698e04984SBarry Smith   PC             cpc;
7573db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
758f3fbd535SBarry Smith   Mat            dA,dB;
759f3fbd535SBarry Smith   Vec            tvec;
760218a07d4SBarry Smith   DM             *dms;
7610a545947SLisandro Dalcin   PetscViewer    viewer = NULL;
76241b6fd38SMatthew G. Knepley   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE;
763f3fbd535SBarry Smith 
764f3fbd535SBarry Smith   PetscFunctionBegin;
7651a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
7661a97d7b6SStefano Zampini   n = mglevels[0]->levels;
76701bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
7683aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
7693aeaf226SBarry Smith     PetscInt levels;
7703aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
7713aeaf226SBarry Smith     levels++;
7723aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
7730298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
7743aeaf226SBarry Smith       n        = levels;
7753aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
7763aeaf226SBarry Smith       mglevels = mg->levels;
7773aeaf226SBarry Smith     }
7783aeaf226SBarry Smith   }
77998e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
7803aeaf226SBarry Smith 
781f3fbd535SBarry Smith 
782f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
783f3fbd535SBarry Smith   /* so use those from global PC */
784f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
7850298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
78698e04984SBarry Smith   if (opsset) {
78798e04984SBarry Smith     Mat mmat;
78823ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
78998e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
79098e04984SBarry Smith   }
7914a5f07a5SMark F. Adams 
79241b6fd38SMatthew G. Knepley   /* Create CR solvers */
79341b6fd38SMatthew G. Knepley   ierr = PCMGGetAdaptCR(pc, &doCR);CHKERRQ(ierr);
79441b6fd38SMatthew G. Knepley   if (doCR) {
79541b6fd38SMatthew G. Knepley     const char *prefix;
79641b6fd38SMatthew G. Knepley 
79741b6fd38SMatthew G. Knepley     ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
79841b6fd38SMatthew G. Knepley     for (i = 1; i < n; ++i) {
79941b6fd38SMatthew G. Knepley       PC   ipc, cr;
80041b6fd38SMatthew G. Knepley       char crprefix[128];
80141b6fd38SMatthew G. Knepley 
80241b6fd38SMatthew G. Knepley       ierr = KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr);CHKERRQ(ierr);
80341b6fd38SMatthew G. Knepley       ierr = KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE);CHKERRQ(ierr);
80441b6fd38SMatthew G. Knepley       ierr = PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i);CHKERRQ(ierr);
80541b6fd38SMatthew G. Knepley       ierr = KSPSetOptionsPrefix(mglevels[i]->cr, prefix);CHKERRQ(ierr);
80641b6fd38SMatthew G. Knepley       ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
80741b6fd38SMatthew G. Knepley       ierr = KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV);CHKERRQ(ierr);
80841b6fd38SMatthew G. Knepley       ierr = KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL);CHKERRQ(ierr);
80941b6fd38SMatthew G. Knepley       ierr = KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED);CHKERRQ(ierr);
81041b6fd38SMatthew G. Knepley       ierr = KSPGetPC(mglevels[i]->cr, &ipc);CHKERRQ(ierr);
81141b6fd38SMatthew G. Knepley 
81241b6fd38SMatthew G. Knepley       ierr = PCSetType(ipc, PCCOMPOSITE);CHKERRQ(ierr);
81341b6fd38SMatthew G. Knepley       ierr = PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
81441b6fd38SMatthew G. Knepley       ierr = PCCompositeAddPCType(ipc, PCSOR);CHKERRQ(ierr);
81541b6fd38SMatthew G. Knepley       ierr = CreateCR_Private(pc, i, &cr);CHKERRQ(ierr);
81641b6fd38SMatthew G. Knepley       ierr = PCCompositeAddPC(ipc, cr);CHKERRQ(ierr);
81741b6fd38SMatthew G. Knepley       ierr = PCDestroy(&cr);CHKERRQ(ierr);
81841b6fd38SMatthew G. Knepley 
81941b6fd38SMatthew G. Knepley       ierr = KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
82041b6fd38SMatthew G. Knepley       ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE);CHKERRQ(ierr);
82141b6fd38SMatthew G. Knepley       ierr = PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i);CHKERRQ(ierr);
82241b6fd38SMatthew G. Knepley       ierr = KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix);CHKERRQ(ierr);
82341b6fd38SMatthew G. Knepley       ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr);CHKERRQ(ierr);
82441b6fd38SMatthew G. Knepley     }
82541b6fd38SMatthew G. Knepley   }
82641b6fd38SMatthew G. Knepley 
8274a5f07a5SMark F. Adams   if (!opsset) {
82871b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
8294a5f07a5SMark F. Adams     if (use_amat) {
830f3fbd535SBarry Smith       ierr = PetscInfo(pc,"Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr);
83123ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
83269858f1bSStefano Zampini     } else {
8334a5f07a5SMark F. Adams       ierr = PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr);
83423ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
8354a5f07a5SMark F. Adams     }
8364a5f07a5SMark F. Adams   }
837f3fbd535SBarry Smith 
8389c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
8399c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
8409c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
8419c7ffb39SBarry Smith       continue;
8429c7ffb39SBarry Smith     }
8439c7ffb39SBarry Smith   }
8442134b1e4SBarry Smith 
8452134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
8462134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
8472134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
8482134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
8492134b1e4SBarry Smith   }
8502134b1e4SBarry Smith 
8512134b1e4SBarry Smith 
8529c7ffb39SBarry Smith   /*
8539c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
8549c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
8559c7ffb39SBarry Smith   */
8562134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
8572d2b81a6SBarry Smith         /* construct the interpolation from the DMs */
8582e499ae9SBarry Smith     Mat p;
85973250ac0SBarry Smith     Vec rscale;
860785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
861218a07d4SBarry Smith     dms[n-1] = pc->dm;
862ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
863ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
86441fce8fdSHong Zhang         /*
86541fce8fdSHong Zhang            Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
86641fce8fdSHong Zhang            Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
86741fce8fdSHong Zhang            But it is safe to use -dm_mat_type.
86841fce8fdSHong Zhang 
86941fce8fdSHong Zhang            The mat type should not be hardcoded like this, we need to find a better way.
87041fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
87141fce8fdSHong Zhang     */
872218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
873942e3340SBarry Smith       DMKSP     kdm;
874eab52d2bSLawrence Mitchell       PetscBool dmhasrestrict, dmhasinject;
8753c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
8762134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
877c27ee7a3SPatrick Farrell       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
878c27ee7a3SPatrick Farrell         ierr = KSPSetDM(mglevels[i]->smoothu,dms[i]);CHKERRQ(ierr);
879c27ee7a3SPatrick Farrell         if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);CHKERRQ(ierr);}
880c27ee7a3SPatrick Farrell       }
88141b6fd38SMatthew G. Knepley       if (mglevels[i]->cr) {
88241b6fd38SMatthew G. Knepley         ierr = KSPSetDM(mglevels[i]->cr,dms[i]);CHKERRQ(ierr);
88341b6fd38SMatthew G. Knepley         if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE);CHKERRQ(ierr);}
88441b6fd38SMatthew G. Knepley       }
885942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
886d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
887d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
8880298fd71SBarry Smith       kdm->ops->computerhs = NULL;
8890298fd71SBarry Smith       kdm->rhsctx          = NULL;
89024c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
891e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
8922d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
89300ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
89473250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
8956bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
896218a07d4SBarry Smith       }
8973ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
8983ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
8993ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
9003ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
9013ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
9023ad4599aSBarry Smith       }
903eab52d2bSLawrence Mitchell       ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr);
904eab52d2bSLawrence Mitchell       if (dmhasinject && !mglevels[i+1]->inject){
905eab52d2bSLawrence Mitchell         ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr);
906eab52d2bSLawrence Mitchell         ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr);
907eab52d2bSLawrence Mitchell         ierr = MatDestroy(&p);CHKERRQ(ierr);
908eab52d2bSLawrence Mitchell       }
90924c3aa18SBarry Smith     }
9102d2b81a6SBarry Smith 
911ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
9122d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
913ce4cda84SJed Brown   }
914cab2e9ccSBarry Smith 
915ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
9162134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
917cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
918cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
919c27ee7a3SPatrick Farrell     if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) {
920c27ee7a3SPatrick Farrell       ierr = KSPSetDM(mglevels[n-1]->smoothu,pc->dm);CHKERRQ(ierr);
921c27ee7a3SPatrick Farrell       ierr = KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);CHKERRQ(ierr);
922c27ee7a3SPatrick Farrell     }
92341b6fd38SMatthew G. Knepley     if (mglevels[n-1]->cr) {
92441b6fd38SMatthew G. Knepley       ierr = KSPSetDM(mglevels[n-1]->cr,pc->dm);CHKERRQ(ierr);
92541b6fd38SMatthew G. Knepley       ierr = KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE);CHKERRQ(ierr);
92641b6fd38SMatthew G. Knepley     }
927218a07d4SBarry Smith   }
928218a07d4SBarry Smith 
9292134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
9302134b1e4SBarry Smith     Mat       A,B;
9312134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
9322134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
9332134b1e4SBarry Smith 
9342134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
9352134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
9362134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
937f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
938b9367914SBarry Smith       if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0");
939b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
940b9367914SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
941b9367914SBarry Smith       }
942b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
943b9367914SBarry Smith         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
944b9367914SBarry Smith       }
9452134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
9462134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
9472134b1e4SBarry Smith       }
9482134b1e4SBarry Smith       if (doA) {
9492df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
9502134b1e4SBarry Smith       }
9512134b1e4SBarry Smith       if (doB) {
9522df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
953b9367914SBarry Smith       }
9542134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
9552134b1e4SBarry Smith       if (!doA && dAeqdB) {
9562134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
9572134b1e4SBarry Smith         A = B;
9582134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
9592134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
9602134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
961b9367914SBarry Smith       }
9622134b1e4SBarry Smith       if (!doB && dAeqdB) {
9632134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
9642134b1e4SBarry Smith         B = A;
9652134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
96623ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
9672134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
9687d537d92SJed Brown       }
9692134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
9702134b1e4SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
9712134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
9722134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
9732134b1e4SBarry Smith       }
9742134b1e4SBarry Smith       dA = A;
975f3fbd535SBarry Smith       dB = B;
976f3fbd535SBarry Smith     }
977f3fbd535SBarry Smith   }
978f3b08a26SMatthew G. Knepley 
979f3b08a26SMatthew G. Knepley 
980f3b08a26SMatthew G. Knepley   /* Adapt interpolation matrices */
981f3b08a26SMatthew G. Knepley   if (mg->adaptInterpolation) {
982f3b08a26SMatthew G. Knepley     mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */
983f3b08a26SMatthew G. Knepley     for (i = 0; i < n; ++i) {
984f3b08a26SMatthew G. Knepley       ierr = PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);CHKERRQ(ierr);
985f3b08a26SMatthew G. Knepley       if (i) {ierr = PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mg->Nc, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace);CHKERRQ(ierr);}
986f3b08a26SMatthew G. Knepley     }
987f3b08a26SMatthew G. Knepley     for (i = n-2; i > -1; --i) {
988f3b08a26SMatthew G. Knepley       ierr = PCMGRecomputeLevelOperators_Internal(pc, i);CHKERRQ(ierr);
989f3b08a26SMatthew G. Knepley     }
990f3b08a26SMatthew G. Knepley   }
991f3b08a26SMatthew G. Knepley 
9922134b1e4SBarry Smith   if (needRestricts && pc->dm) {
993caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
994ccceb50cSJed Brown       DM  dmfine,dmcoarse;
995caa4e7f2SJed Brown       Mat Restrict,Inject;
996caa4e7f2SJed Brown       Vec rscale;
997ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
998ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
999caa4e7f2SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
1000caa4e7f2SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
1001eab52d2bSLawrence Mitchell       ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr);
1002caa4e7f2SJed Brown       ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
1003caa4e7f2SJed Brown     }
1004caa4e7f2SJed Brown   }
1005f3fbd535SBarry Smith 
1006f3fbd535SBarry Smith   if (!pc->setupcalled) {
1007f3fbd535SBarry Smith     for (i=0; i<n; i++) {
1008f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
1009f3fbd535SBarry Smith     }
1010f3fbd535SBarry Smith     for (i=1; i<n; i++) {
1011f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
1012f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
1013f3fbd535SBarry Smith       }
101441b6fd38SMatthew G. Knepley       if (mglevels[i]->cr) {
101541b6fd38SMatthew G. Knepley         ierr = KSPSetFromOptions(mglevels[i]->cr);CHKERRQ(ierr);
101641b6fd38SMatthew G. Knepley       }
1017f3fbd535SBarry Smith     }
10183ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
1019f3fbd535SBarry Smith     for (i=1; i<n; i++) {
10203ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
10213ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
1022f3fbd535SBarry Smith     }
1023f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
1024f3fbd535SBarry Smith       if (!mglevels[i]->b) {
1025f3fbd535SBarry Smith         Vec *vec;
10262a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
1027f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
10286bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
1029f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
1030f3fbd535SBarry Smith       }
1031f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
1032f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
1033f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
10346bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
1035f3fbd535SBarry Smith       }
1036f3fbd535SBarry Smith       if (!mglevels[i]->x) {
1037f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
1038f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
10396bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
1040f3fbd535SBarry Smith       }
104141b6fd38SMatthew G. Knepley       if (doCR) {
104241b6fd38SMatthew G. Knepley         ierr = VecDuplicate(mglevels[i]->b,&mglevels[i]->crx);CHKERRQ(ierr);
104341b6fd38SMatthew G. Knepley         ierr = VecDuplicate(mglevels[i]->b,&mglevels[i]->crb);CHKERRQ(ierr);
104441b6fd38SMatthew G. Knepley         ierr = VecZeroEntries(mglevels[i]->crb);CHKERRQ(ierr);
104541b6fd38SMatthew G. Knepley       }
1046f3fbd535SBarry Smith     }
1047f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
1048f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
1049f3fbd535SBarry Smith       Vec *vec;
10502a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
1051f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
10526bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
1053f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
1054f3fbd535SBarry Smith     }
105541b6fd38SMatthew G. Knepley     if (doCR) {
105641b6fd38SMatthew G. Knepley       ierr = VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx);CHKERRQ(ierr);
105741b6fd38SMatthew G. Knepley       ierr = VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb);CHKERRQ(ierr);
105841b6fd38SMatthew G. Knepley       ierr = VecZeroEntries(mglevels[n-1]->crb);CHKERRQ(ierr);
105941b6fd38SMatthew G. Knepley     }
1060f3fbd535SBarry Smith   }
1061f3fbd535SBarry Smith 
106298e04984SBarry Smith   if (pc->dm) {
106398e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
106498e04984SBarry Smith     for (i=0; i<n-1; i++) {
106598e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
106698e04984SBarry Smith     }
106798e04984SBarry Smith   }
1068f3fbd535SBarry Smith 
1069f3fbd535SBarry Smith   for (i=1; i<n; i++) {
10702a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
1071f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
1072f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
1073f3fbd535SBarry Smith     }
107441b6fd38SMatthew G. Knepley     if (mglevels[i]->cr) {ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);CHKERRQ(ierr);}
107563e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1076f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
1077c0decd05SBarry Smith     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1078899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1079899639b0SHong Zhang     }
108063e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1081d42688cbSBarry Smith     if (!mglevels[i]->residual) {
1082d42688cbSBarry Smith       Mat mat;
108313842ffbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
108454b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
1085d42688cbSBarry Smith     }
1086f3fbd535SBarry Smith   }
1087f3fbd535SBarry Smith   for (i=1; i<n; i++) {
1088f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
1089f3fbd535SBarry Smith       Mat downmat,downpmat;
1090f3fbd535SBarry Smith 
1091f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
10920298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
1093f3fbd535SBarry Smith       if (!opsset) {
109423ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
109523ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
1096f3fbd535SBarry Smith       }
1097f3fbd535SBarry Smith 
1098f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
109963e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1100f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
1101c0decd05SBarry Smith       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
1102899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
1103899639b0SHong Zhang       }
110463e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1105f3fbd535SBarry Smith     }
110641b6fd38SMatthew G. Knepley     if (mglevels[i]->cr) {
110741b6fd38SMatthew G. Knepley       Mat downmat,downpmat;
110841b6fd38SMatthew G. Knepley 
110941b6fd38SMatthew G. Knepley       /* check if operators have been set for up, if not use down operators to set them */
111041b6fd38SMatthew G. Knepley       ierr = KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL);CHKERRQ(ierr);
111141b6fd38SMatthew G. Knepley       if (!opsset) {
111241b6fd38SMatthew G. Knepley         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
111341b6fd38SMatthew G. Knepley         ierr = KSPSetOperators(mglevels[i]->cr,downmat,downpmat);CHKERRQ(ierr);
111441b6fd38SMatthew G. Knepley       }
111541b6fd38SMatthew G. Knepley 
111641b6fd38SMatthew G. Knepley       ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);CHKERRQ(ierr);
111741b6fd38SMatthew G. Knepley       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
111841b6fd38SMatthew G. Knepley       ierr = KSPSetUp(mglevels[i]->cr);CHKERRQ(ierr);
111941b6fd38SMatthew G. Knepley       if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) {
112041b6fd38SMatthew G. Knepley         pc->failedreason = PC_SUBPC_ERROR;
112141b6fd38SMatthew G. Knepley       }
112241b6fd38SMatthew G. Knepley       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
112341b6fd38SMatthew G. Knepley     }
1124f3fbd535SBarry Smith   }
1125f3fbd535SBarry Smith 
112663e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1127f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
1128c0decd05SBarry Smith   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1129899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
1130899639b0SHong Zhang   }
113163e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1132f3fbd535SBarry Smith 
1133f3fbd535SBarry Smith   /*
1134f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
1135e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
1136f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
1137f3fbd535SBarry Smith 
1138f3fbd535SBarry Smith    Only support one or the other at the same time.
1139f3fbd535SBarry Smith   */
1140f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
1141c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
1142ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
1143f3fbd535SBarry Smith   dump = PETSC_FALSE;
1144f3fbd535SBarry Smith #endif
1145c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
1146ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
1147f3fbd535SBarry Smith 
1148f3fbd535SBarry Smith   if (viewer) {
1149f3fbd535SBarry Smith     for (i=1; i<n; i++) {
1150f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
1151f3fbd535SBarry Smith     }
1152f3fbd535SBarry Smith     for (i=0; i<n; i++) {
1153f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
1154f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1155f3fbd535SBarry Smith     }
1156f3fbd535SBarry Smith   }
1157f3fbd535SBarry Smith   PetscFunctionReturn(0);
1158f3fbd535SBarry Smith }
1159f3fbd535SBarry Smith 
1160f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
1161f3fbd535SBarry Smith 
116297d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
116397d33e41SMatthew G. Knepley {
116497d33e41SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
116597d33e41SMatthew G. Knepley 
116697d33e41SMatthew G. Knepley   PetscFunctionBegin;
116797d33e41SMatthew G. Knepley   *levels = mg->nlevels;
116897d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
116997d33e41SMatthew G. Knepley }
117097d33e41SMatthew G. Knepley 
11714b9ad928SBarry Smith /*@
117297177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
11734b9ad928SBarry Smith 
11744b9ad928SBarry Smith    Not Collective
11754b9ad928SBarry Smith 
11764b9ad928SBarry Smith    Input Parameter:
11774b9ad928SBarry Smith .  pc - the preconditioner context
11784b9ad928SBarry Smith 
11794b9ad928SBarry Smith    Output parameter:
11804b9ad928SBarry Smith .  levels - the number of levels
11814b9ad928SBarry Smith 
11824b9ad928SBarry Smith    Level: advanced
11834b9ad928SBarry Smith 
118497177400SBarry Smith .seealso: PCMGSetLevels()
11854b9ad928SBarry Smith @*/
11867087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
11874b9ad928SBarry Smith {
1188e952c529SMatthew G. Knepley   PetscErrorCode ierr;
11894b9ad928SBarry Smith 
11904b9ad928SBarry Smith   PetscFunctionBegin;
11910700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11924482741eSBarry Smith   PetscValidIntPointer(levels,2);
119397d33e41SMatthew G. Knepley   *levels = 0;
119437b5128cSMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr);
11954b9ad928SBarry Smith   PetscFunctionReturn(0);
11964b9ad928SBarry Smith }
11974b9ad928SBarry Smith 
11984b9ad928SBarry Smith /*@
119997177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
12004b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
12014b9ad928SBarry Smith 
1202ad4df100SBarry Smith    Logically Collective on PC
12034b9ad928SBarry Smith 
12044b9ad928SBarry Smith    Input Parameters:
12054b9ad928SBarry Smith +  pc - the preconditioner context
12069dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
12079dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
12084b9ad928SBarry Smith 
12094b9ad928SBarry Smith    Options Database Key:
12104b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
12114b9ad928SBarry Smith    additive, full, kaskade
12124b9ad928SBarry Smith 
12134b9ad928SBarry Smith    Level: advanced
12144b9ad928SBarry Smith 
121597177400SBarry Smith .seealso: PCMGSetLevels()
12164b9ad928SBarry Smith @*/
12177087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
12184b9ad928SBarry Smith {
1219f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
12204b9ad928SBarry Smith 
12214b9ad928SBarry Smith   PetscFunctionBegin;
12220700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1223c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
122431567311SBarry Smith   mg->am = form;
12259dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1226c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
1227c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1228c60c7ad4SBarry Smith }
1229c60c7ad4SBarry Smith 
1230c60c7ad4SBarry Smith /*@
1231c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
1232c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
1233c60c7ad4SBarry Smith 
1234c60c7ad4SBarry Smith    Logically Collective on PC
1235c60c7ad4SBarry Smith 
1236c60c7ad4SBarry Smith    Input Parameter:
1237c60c7ad4SBarry Smith .  pc - the preconditioner context
1238c60c7ad4SBarry Smith 
1239c60c7ad4SBarry Smith    Output Parameter:
1240c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
1241c60c7ad4SBarry Smith 
1242c60c7ad4SBarry Smith 
1243c60c7ad4SBarry Smith    Level: advanced
1244c60c7ad4SBarry Smith 
1245c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
1246c60c7ad4SBarry Smith @*/
1247c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1248c60c7ad4SBarry Smith {
1249c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1250c60c7ad4SBarry Smith 
1251c60c7ad4SBarry Smith   PetscFunctionBegin;
1252c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1253c60c7ad4SBarry Smith   *type = mg->am;
12544b9ad928SBarry Smith   PetscFunctionReturn(0);
12554b9ad928SBarry Smith }
12564b9ad928SBarry Smith 
12574b9ad928SBarry Smith /*@
12580d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
12594b9ad928SBarry Smith    complicated cycling.
12604b9ad928SBarry Smith 
1261ad4df100SBarry Smith    Logically Collective on PC
12624b9ad928SBarry Smith 
12634b9ad928SBarry Smith    Input Parameters:
1264c2be2410SBarry Smith +  pc - the multigrid context
1265c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
12664b9ad928SBarry Smith 
12674b9ad928SBarry Smith    Options Database Key:
1268c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
12694b9ad928SBarry Smith 
12704b9ad928SBarry Smith    Level: advanced
12714b9ad928SBarry Smith 
12720d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
12734b9ad928SBarry Smith @*/
12747087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
12754b9ad928SBarry Smith {
1276f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
1277f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
127879416396SBarry Smith   PetscInt     i,levels;
12794b9ad928SBarry Smith 
12804b9ad928SBarry Smith   PetscFunctionBegin;
12810700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
128269fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
12831a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1284f3fbd535SBarry Smith   levels = mglevels[0]->levels;
12852fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
12864b9ad928SBarry Smith   PetscFunctionReturn(0);
12874b9ad928SBarry Smith }
12884b9ad928SBarry Smith 
12898cc2d5dfSBarry Smith /*@
12908cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
12918cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
12928cc2d5dfSBarry Smith 
1293ad4df100SBarry Smith    Logically Collective on PC
12948cc2d5dfSBarry Smith 
12958cc2d5dfSBarry Smith    Input Parameters:
12968cc2d5dfSBarry Smith +  pc - the multigrid context
12978cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
12988cc2d5dfSBarry Smith 
12998cc2d5dfSBarry Smith    Options Database Key:
1300e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
13018cc2d5dfSBarry Smith 
13028cc2d5dfSBarry Smith    Level: advanced
13038cc2d5dfSBarry Smith 
130495452b02SPatrick Sanan    Notes:
130595452b02SPatrick Sanan     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
13068cc2d5dfSBarry Smith 
13078cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
13088cc2d5dfSBarry Smith @*/
13097087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
13108cc2d5dfSBarry Smith {
1311f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
13128cc2d5dfSBarry Smith 
13138cc2d5dfSBarry Smith   PetscFunctionBegin;
13140700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1315c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
13162a21e185SBarry Smith   mg->cyclesperpcapply = n;
13178cc2d5dfSBarry Smith   PetscFunctionReturn(0);
13188cc2d5dfSBarry Smith }
13198cc2d5dfSBarry Smith 
13202134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1321fb15c04eSBarry Smith {
1322fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1323fb15c04eSBarry Smith 
1324fb15c04eSBarry Smith   PetscFunctionBegin;
13252134b1e4SBarry Smith   mg->galerkin = use;
1326fb15c04eSBarry Smith   PetscFunctionReturn(0);
1327fb15c04eSBarry Smith }
1328fb15c04eSBarry Smith 
1329c2be2410SBarry Smith /*@
133097177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
133182b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1332c2be2410SBarry Smith 
1333ad4df100SBarry Smith    Logically Collective on PC
1334c2be2410SBarry Smith 
1335c2be2410SBarry Smith    Input Parameters:
1336c91913d3SJed Brown +  pc - the multigrid context
13372134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1338c2be2410SBarry Smith 
1339c2be2410SBarry Smith    Options Database Key:
13402134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1341c2be2410SBarry Smith 
1342c2be2410SBarry Smith    Level: intermediate
1343c2be2410SBarry Smith 
134495452b02SPatrick Sanan    Notes:
134595452b02SPatrick Sanan     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
13461c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
13471c1aac46SBarry Smith 
13482134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
13493fc8bf9cSBarry Smith 
1350c2be2410SBarry Smith @*/
13512134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1352c2be2410SBarry Smith {
1353fb15c04eSBarry Smith   PetscErrorCode ierr;
1354c2be2410SBarry Smith 
1355c2be2410SBarry Smith   PetscFunctionBegin;
13560700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13572134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1358c2be2410SBarry Smith   PetscFunctionReturn(0);
1359c2be2410SBarry Smith }
1360c2be2410SBarry Smith 
13613fc8bf9cSBarry Smith /*@
13623fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
136382b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
13643fc8bf9cSBarry Smith 
13653fc8bf9cSBarry Smith    Not Collective
13663fc8bf9cSBarry Smith 
13673fc8bf9cSBarry Smith    Input Parameter:
13683fc8bf9cSBarry Smith .  pc - the multigrid context
13693fc8bf9cSBarry Smith 
13703fc8bf9cSBarry Smith    Output Parameter:
13712134b1e4SBarry 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
13723fc8bf9cSBarry Smith 
13733fc8bf9cSBarry Smith    Level: intermediate
13743fc8bf9cSBarry Smith 
13752134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
13763fc8bf9cSBarry Smith 
13773fc8bf9cSBarry Smith @*/
13782134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
13793fc8bf9cSBarry Smith {
1380f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
13813fc8bf9cSBarry Smith 
13823fc8bf9cSBarry Smith   PetscFunctionBegin;
13830700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13842134b1e4SBarry Smith   *galerkin = mg->galerkin;
13853fc8bf9cSBarry Smith   PetscFunctionReturn(0);
13863fc8bf9cSBarry Smith }
13873fc8bf9cSBarry Smith 
1388f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1389f3b08a26SMatthew G. Knepley {
1390f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
1391f3b08a26SMatthew G. Knepley 
1392f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1393f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = adapt;
1394f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1395f3b08a26SMatthew G. Knepley }
1396f3b08a26SMatthew G. Knepley 
1397f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1398f3b08a26SMatthew G. Knepley {
1399f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
1400f3b08a26SMatthew G. Knepley 
1401f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1402f3b08a26SMatthew G. Knepley   *adapt = mg->adaptInterpolation;
1403f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1404f3b08a26SMatthew G. Knepley }
1405f3b08a26SMatthew G. Knepley 
140641b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr)
140741b6fd38SMatthew G. Knepley {
140841b6fd38SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
140941b6fd38SMatthew G. Knepley 
141041b6fd38SMatthew G. Knepley   PetscFunctionBegin;
141141b6fd38SMatthew G. Knepley   mg->compatibleRelaxation = cr;
141241b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
141341b6fd38SMatthew G. Knepley }
141441b6fd38SMatthew G. Knepley 
141541b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr)
141641b6fd38SMatthew G. Knepley {
141741b6fd38SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
141841b6fd38SMatthew G. Knepley 
141941b6fd38SMatthew G. Knepley   PetscFunctionBegin;
142041b6fd38SMatthew G. Knepley   *cr = mg->compatibleRelaxation;
142141b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
142241b6fd38SMatthew G. Knepley }
142341b6fd38SMatthew G. Knepley 
1424f3b08a26SMatthew G. Knepley /*@
1425f3b08a26SMatthew 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.
1426f3b08a26SMatthew G. Knepley 
1427f3b08a26SMatthew G. Knepley   Logically Collective on PC
1428f3b08a26SMatthew G. Knepley 
1429f3b08a26SMatthew G. Knepley   Input Parameters:
1430f3b08a26SMatthew G. Knepley + pc    - the multigrid context
1431f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator
1432f3b08a26SMatthew G. Knepley 
1433f3b08a26SMatthew G. Knepley   Options Database Keys:
1434f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp                     - Turn on adaptation
1435f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int>             - The number of modes to use, should be divisible by dimension
1436f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector
1437f3b08a26SMatthew G. Knepley 
1438f3b08a26SMatthew G. Knepley   Level: intermediate
1439f3b08a26SMatthew G. Knepley 
1440f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin
1441f3b08a26SMatthew G. Knepley .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1442f3b08a26SMatthew G. Knepley @*/
1443f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1444f3b08a26SMatthew G. Knepley {
1445f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1446f3b08a26SMatthew G. Knepley 
1447f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1448f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1449f3b08a26SMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));CHKERRQ(ierr);
1450f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1451f3b08a26SMatthew G. Knepley }
1452f3b08a26SMatthew G. Knepley 
1453f3b08a26SMatthew G. Knepley /*@
1454f3b08a26SMatthew 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.
1455f3b08a26SMatthew G. Knepley 
1456f3b08a26SMatthew G. Knepley   Not collective
1457f3b08a26SMatthew G. Knepley 
1458f3b08a26SMatthew G. Knepley   Input Parameter:
1459f3b08a26SMatthew G. Knepley . pc    - the multigrid context
1460f3b08a26SMatthew G. Knepley 
1461f3b08a26SMatthew G. Knepley   Output Parameter:
1462f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator
1463f3b08a26SMatthew G. Knepley 
1464f3b08a26SMatthew G. Knepley   Level: intermediate
1465f3b08a26SMatthew G. Knepley 
1466f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin
1467f3b08a26SMatthew G. Knepley .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1468f3b08a26SMatthew G. Knepley @*/
1469f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1470f3b08a26SMatthew G. Knepley {
1471f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1472f3b08a26SMatthew G. Knepley 
1473f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1474f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1475f3b08a26SMatthew G. Knepley   PetscValidPointer(adapt, 2);
1476f3b08a26SMatthew G. Knepley   ierr = PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));CHKERRQ(ierr);
1477f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1478f3b08a26SMatthew G. Knepley }
1479f3b08a26SMatthew G. Knepley 
14804b9ad928SBarry Smith /*@
148141b6fd38SMatthew G. Knepley   PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation.
148241b6fd38SMatthew G. Knepley 
148341b6fd38SMatthew G. Knepley   Logically Collective on PC
148441b6fd38SMatthew G. Knepley 
148541b6fd38SMatthew G. Knepley   Input Parameters:
148641b6fd38SMatthew G. Knepley + pc - the multigrid context
148741b6fd38SMatthew G. Knepley - cr - flag for compatible relaxation
148841b6fd38SMatthew G. Knepley 
148941b6fd38SMatthew G. Knepley   Options Database Keys:
149041b6fd38SMatthew G. Knepley . -pc_mg_adapt_cr - Turn on compatible relaxation
149141b6fd38SMatthew G. Knepley 
149241b6fd38SMatthew G. Knepley   Level: intermediate
149341b6fd38SMatthew G. Knepley 
149441b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin
149541b6fd38SMatthew G. Knepley .seealso: PCMGGetAdaptCR(), PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
149641b6fd38SMatthew G. Knepley @*/
149741b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr)
149841b6fd38SMatthew G. Knepley {
149941b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
150041b6fd38SMatthew G. Knepley 
150141b6fd38SMatthew G. Knepley   PetscFunctionBegin;
150241b6fd38SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
150341b6fd38SMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr));CHKERRQ(ierr);
150441b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
150541b6fd38SMatthew G. Knepley }
150641b6fd38SMatthew G. Knepley 
150741b6fd38SMatthew G. Knepley /*@
150841b6fd38SMatthew G. Knepley   PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation.
150941b6fd38SMatthew G. Knepley 
151041b6fd38SMatthew G. Knepley   Not collective
151141b6fd38SMatthew G. Knepley 
151241b6fd38SMatthew G. Knepley   Input Parameter:
151341b6fd38SMatthew G. Knepley . pc    - the multigrid context
151441b6fd38SMatthew G. Knepley 
151541b6fd38SMatthew G. Knepley   Output Parameter:
151641b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion
151741b6fd38SMatthew G. Knepley 
151841b6fd38SMatthew G. Knepley   Level: intermediate
151941b6fd38SMatthew G. Knepley 
152041b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin
152141b6fd38SMatthew G. Knepley .seealso: PCMGSetAdaptCR(), PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
152241b6fd38SMatthew G. Knepley @*/
152341b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr)
152441b6fd38SMatthew G. Knepley {
152541b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
152641b6fd38SMatthew G. Knepley 
152741b6fd38SMatthew G. Knepley   PetscFunctionBegin;
152841b6fd38SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
152941b6fd38SMatthew G. Knepley   PetscValidPointer(cr, 2);
153041b6fd38SMatthew G. Knepley   ierr = PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr));CHKERRQ(ierr);
153141b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
153241b6fd38SMatthew G. Knepley }
153341b6fd38SMatthew G. Knepley 
153441b6fd38SMatthew G. Knepley /*@
153506792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1536710315b6SLawrence Mitchell    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1537710315b6SLawrence Mitchell    pre- and post-smoothing steps.
153806792cafSBarry Smith 
153906792cafSBarry Smith    Logically Collective on PC
154006792cafSBarry Smith 
154106792cafSBarry Smith    Input Parameters:
154206792cafSBarry Smith +  mg - the multigrid context
154306792cafSBarry Smith -  n - the number of smoothing steps
154406792cafSBarry Smith 
154506792cafSBarry Smith    Options Database Key:
1546a2b725a8SWilliam Gropp .  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
154706792cafSBarry Smith 
154806792cafSBarry Smith    Level: advanced
154906792cafSBarry Smith 
155095452b02SPatrick Sanan    Notes:
155195452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
155206792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
155306792cafSBarry Smith 
1554710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp()
155506792cafSBarry Smith @*/
155606792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
155706792cafSBarry Smith {
155806792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
155906792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
156006792cafSBarry Smith   PetscErrorCode ierr;
156106792cafSBarry Smith   PetscInt       i,levels;
156206792cafSBarry Smith 
156306792cafSBarry Smith   PetscFunctionBegin;
156406792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
156506792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
15661a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
156706792cafSBarry Smith   levels = mglevels[0]->levels;
156806792cafSBarry Smith 
156906792cafSBarry Smith   for (i=1; i<levels; i++) {
157006792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
157106792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
157206792cafSBarry Smith     mg->default_smoothu = n;
157306792cafSBarry Smith     mg->default_smoothd = n;
157406792cafSBarry Smith   }
157506792cafSBarry Smith   PetscFunctionReturn(0);
157606792cafSBarry Smith }
157706792cafSBarry Smith 
1578f442ab6aSBarry Smith /*@
15798e5aa403SBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels
1580710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1581f442ab6aSBarry Smith 
1582f442ab6aSBarry Smith    Logically Collective on PC
1583f442ab6aSBarry Smith 
1584f442ab6aSBarry Smith    Input Parameters:
1585f442ab6aSBarry Smith .  pc - the preconditioner context
1586f442ab6aSBarry Smith 
1587f442ab6aSBarry Smith    Options Database Key:
1588f442ab6aSBarry Smith .  -pc_mg_distinct_smoothup
1589f442ab6aSBarry Smith 
1590f442ab6aSBarry Smith    Level: advanced
1591f442ab6aSBarry Smith 
159295452b02SPatrick Sanan    Notes:
159395452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
1594f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1595f442ab6aSBarry Smith 
1596710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth()
1597f442ab6aSBarry Smith @*/
1598f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1599f442ab6aSBarry Smith {
1600f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1601f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1602f442ab6aSBarry Smith   PetscErrorCode ierr;
1603f442ab6aSBarry Smith   PetscInt       i,levels;
1604f442ab6aSBarry Smith   KSP            subksp;
1605f442ab6aSBarry Smith 
1606f442ab6aSBarry Smith   PetscFunctionBegin;
1607f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1608f442ab6aSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1609f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1610f442ab6aSBarry Smith 
1611f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1612710315b6SLawrence Mitchell     const char *prefix = NULL;
1613f442ab6aSBarry Smith     /* make sure smoother up and down are different */
1614f442ab6aSBarry Smith     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1615710315b6SLawrence Mitchell     ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr);
1616710315b6SLawrence Mitchell     ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr);
1617f442ab6aSBarry Smith     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1618f442ab6aSBarry Smith   }
1619f442ab6aSBarry Smith   PetscFunctionReturn(0);
1620f442ab6aSBarry Smith }
1621f442ab6aSBarry Smith 
162207a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
162307a4832bSFande Kong PetscErrorCode  PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[])
162407a4832bSFande Kong {
162507a4832bSFande Kong   PC_MG          *mg        = (PC_MG*)pc->data;
162607a4832bSFande Kong   PC_MG_Levels   **mglevels = mg->levels;
162707a4832bSFande Kong   Mat            *mat;
162807a4832bSFande Kong   PetscInt       l;
162907a4832bSFande Kong   PetscErrorCode ierr;
163007a4832bSFande Kong 
163107a4832bSFande Kong   PetscFunctionBegin;
163207a4832bSFande Kong   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
163307a4832bSFande Kong   ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr);
163407a4832bSFande Kong   for (l=1; l< mg->nlevels; l++) {
163507a4832bSFande Kong     mat[l-1] = mglevels[l]->interpolate;
163607a4832bSFande Kong     ierr = PetscObjectReference((PetscObject)mat[l-1]);CHKERRQ(ierr);
163707a4832bSFande Kong   }
163807a4832bSFande Kong   *num_levels = mg->nlevels;
163907a4832bSFande Kong   *interpolations = mat;
164007a4832bSFande Kong   PetscFunctionReturn(0);
164107a4832bSFande Kong }
164207a4832bSFande Kong 
164307a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
164407a4832bSFande Kong PetscErrorCode  PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[])
164507a4832bSFande Kong {
164607a4832bSFande Kong   PC_MG          *mg        = (PC_MG*)pc->data;
164707a4832bSFande Kong   PC_MG_Levels   **mglevels = mg->levels;
164807a4832bSFande Kong   PetscInt       l;
164907a4832bSFande Kong   Mat            *mat;
165007a4832bSFande Kong   PetscErrorCode ierr;
165107a4832bSFande Kong 
165207a4832bSFande Kong   PetscFunctionBegin;
165307a4832bSFande Kong   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
165407a4832bSFande Kong   ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr);
165507a4832bSFande Kong   for (l=0; l<mg->nlevels-1; l++) {
165607a4832bSFande Kong     ierr = KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));CHKERRQ(ierr);
165707a4832bSFande Kong     ierr = PetscObjectReference((PetscObject)mat[l]);CHKERRQ(ierr);
165807a4832bSFande Kong   }
165907a4832bSFande Kong   *num_levels = mg->nlevels;
166007a4832bSFande Kong   *coarseOperators = mat;
166107a4832bSFande Kong   PetscFunctionReturn(0);
166207a4832bSFande Kong }
166307a4832bSFande Kong 
1664f3b08a26SMatthew G. Knepley /*@C
1665f3b08a26SMatthew G. Knepley   PCMGRegisterCoarseSpaceConstructor -  Adds a method to the PCMG package for coarse space construction.
1666f3b08a26SMatthew G. Knepley 
1667f3b08a26SMatthew G. Knepley   Not collective
1668f3b08a26SMatthew G. Knepley 
1669f3b08a26SMatthew G. Knepley   Input Parameters:
1670f3b08a26SMatthew G. Knepley + name     - name of the constructor
1671f3b08a26SMatthew G. Knepley - function - constructor routine
1672f3b08a26SMatthew G. Knepley 
1673f3b08a26SMatthew G. Knepley   Notes:
1674f3b08a26SMatthew G. Knepley   Calling sequence for the routine:
1675f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1676f3b08a26SMatthew G. Knepley $   pc        - The PC object
1677f3b08a26SMatthew G. Knepley $   l         - The multigrid level, 0 is the coarse level
1678f3b08a26SMatthew G. Knepley $   dm        - The DM for this level
1679f3b08a26SMatthew G. Knepley $   smooth    - The level smoother
1680f3b08a26SMatthew G. Knepley $   Nc        - The size of the coarse space
1681f3b08a26SMatthew G. Knepley $   initGuess - Basis for an initial guess for the space
1682f3b08a26SMatthew G. Knepley $   coarseSp  - A basis for the computed coarse space
1683f3b08a26SMatthew G. Knepley 
1684f3b08a26SMatthew G. Knepley   Level: advanced
1685f3b08a26SMatthew G. Knepley 
1686f3b08a26SMatthew G. Knepley .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister()
1687f3b08a26SMatthew G. Knepley @*/
1688f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1689f3b08a26SMatthew G. Knepley {
1690f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1691f3b08a26SMatthew G. Knepley 
1692f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1693f3b08a26SMatthew G. Knepley   ierr = PCInitializePackage();CHKERRQ(ierr);
1694f3b08a26SMatthew G. Knepley   ierr = PetscFunctionListAdd(&PCMGCoarseList,name,function);CHKERRQ(ierr);
1695f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1696f3b08a26SMatthew G. Knepley }
1697f3b08a26SMatthew G. Knepley 
1698f3b08a26SMatthew G. Knepley /*@C
1699f3b08a26SMatthew G. Knepley   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.
1700f3b08a26SMatthew G. Knepley 
1701f3b08a26SMatthew G. Knepley   Not collective
1702f3b08a26SMatthew G. Knepley 
1703f3b08a26SMatthew G. Knepley   Input Parameter:
1704f3b08a26SMatthew G. Knepley . name     - name of the constructor
1705f3b08a26SMatthew G. Knepley 
1706f3b08a26SMatthew G. Knepley   Output Parameter:
1707f3b08a26SMatthew G. Knepley . function - constructor routine
1708f3b08a26SMatthew G. Knepley 
1709f3b08a26SMatthew G. Knepley   Notes:
1710f3b08a26SMatthew G. Knepley   Calling sequence for the routine:
1711f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1712f3b08a26SMatthew G. Knepley $   pc        - The PC object
1713f3b08a26SMatthew G. Knepley $   l         - The multigrid level, 0 is the coarse level
1714f3b08a26SMatthew G. Knepley $   dm        - The DM for this level
1715f3b08a26SMatthew G. Knepley $   smooth    - The level smoother
1716f3b08a26SMatthew G. Knepley $   Nc        - The size of the coarse space
1717f3b08a26SMatthew G. Knepley $   initGuess - Basis for an initial guess for the space
1718f3b08a26SMatthew G. Knepley $   coarseSp  - A basis for the computed coarse space
1719f3b08a26SMatthew G. Knepley 
1720f3b08a26SMatthew G. Knepley   Level: advanced
1721f3b08a26SMatthew G. Knepley 
1722f3b08a26SMatthew G. Knepley .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister()
1723f3b08a26SMatthew G. Knepley @*/
1724f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1725f3b08a26SMatthew G. Knepley {
1726f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1727f3b08a26SMatthew G. Knepley 
1728f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1729f3b08a26SMatthew G. Knepley   ierr = PetscFunctionListFind(PCMGCoarseList,name,function);CHKERRQ(ierr);
1730f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1731f3b08a26SMatthew G. Knepley }
1732f3b08a26SMatthew G. Knepley 
17334b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
17344b9ad928SBarry Smith 
17353b09bd56SBarry Smith /*MC
1736ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
17373b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
17383b09bd56SBarry Smith 
17393b09bd56SBarry Smith    Options Database Keys:
17403b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1741391689abSStefano Zampini .  -pc_mg_cycle_type <v,w> - provide the cycle desired
17428c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
17433b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1744710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
17452134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
17468cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
17478cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1748e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
17498cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
17508cf6037eSBarry Smith                         to the binary output file called binaryoutput
17513b09bd56SBarry Smith 
175295452b02SPatrick Sanan    Notes:
17530b3c753eSRichard 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
17543b09bd56SBarry Smith 
17558cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
17568cf6037eSBarry Smith 
175723067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
175823067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
175923067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
176023067569SBarry Smith        residual is computed at the end of each cycle.
176123067569SBarry Smith 
17623b09bd56SBarry Smith    Level: intermediate
17633b09bd56SBarry Smith 
17648cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1765710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1766710315b6SLawrence Mitchell            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
176797177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
17680d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
17693b09bd56SBarry Smith M*/
17703b09bd56SBarry Smith 
17718cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
17724b9ad928SBarry Smith {
1773f3fbd535SBarry Smith   PC_MG          *mg;
1774f3fbd535SBarry Smith   PetscErrorCode ierr;
1775f3fbd535SBarry Smith 
17764b9ad928SBarry Smith   PetscFunctionBegin;
1777b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1778f3fbd535SBarry Smith   pc->data     = (void*)mg;
1779f3fbd535SBarry Smith   mg->nlevels  = -1;
178010eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
17812134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1782f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = PETSC_FALSE;
1783f3b08a26SMatthew G. Knepley   mg->Nc                 = -1;
1784f3b08a26SMatthew G. Knepley   mg->eigenvalue         = -1;
1785f3fbd535SBarry Smith 
178637a44384SMark Adams   pc->useAmat = PETSC_TRUE;
178737a44384SMark Adams 
17884b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
17894b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1790a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
17914b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
17924b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
17934b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1794fb15c04eSBarry Smith 
1795f3b08a26SMatthew G. Knepley   ierr = PetscObjectComposedDataRegister(&mg->eigenvalue);CHKERRQ(ierr);
1796fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
179797d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
179897d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
1799fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);CHKERRQ(ierr);
1800fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);CHKERRQ(ierr);
1801f3b08a26SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);CHKERRQ(ierr);
1802f3b08a26SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);CHKERRQ(ierr);
180341b6fd38SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG);CHKERRQ(ierr);
180441b6fd38SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG);CHKERRQ(ierr);
18054b9ad928SBarry Smith   PetscFunctionReturn(0);
18064b9ad928SBarry Smith }
1807