xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision fcb023d4cd21af6e7bfc4dfb4b9546d75029f707)
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 
15*fcb023d4SJed Brown PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PetscBool transpose,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);}
24*fcb023d4SJed Brown   if (!transpose) {
2531567311SBarry Smith     ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
26c0decd05SBarry Smith     ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr);
27*fcb023d4SJed Brown   } else {
28*fcb023d4SJed Brown     ierr = KSPSolveTranspose(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* transpose of post-smooth */
29*fcb023d4SJed Brown     ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr);
30*fcb023d4SJed Brown   }
3163e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
3231567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
3363e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
34*fcb023d4SJed Brown     if (!transpose) {
3531567311SBarry Smith       ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
36*fcb023d4SJed Brown     } else {
37*fcb023d4SJed Brown       ierr = (*mglevels->residualtranspose)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
38*fcb023d4SJed Brown     }
3963e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
404b9ad928SBarry Smith 
414b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
4231567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
434b9ad928SBarry Smith       PetscReal rnorm;
4431567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
454b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
4670441072SBarry Smith         if (rnorm < mg->abstol) {
474d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
4857622a8eSBarry 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);
494b9ad928SBarry Smith         } else {
504d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
5157622a8eSBarry 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);
524b9ad928SBarry Smith         }
534b9ad928SBarry Smith         PetscFunctionReturn(0);
544b9ad928SBarry Smith       }
554b9ad928SBarry Smith     }
564b9ad928SBarry Smith 
5731567311SBarry Smith     mgc = *(mglevelsin - 1);
5863e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
59*fcb023d4SJed Brown     if (!transpose) {
6031567311SBarry Smith       ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
61*fcb023d4SJed Brown     } else {
62*fcb023d4SJed Brown       ierr = MatRestrict(mglevels->interpolate,mglevels->r,mgc->b);CHKERRQ(ierr);
63*fcb023d4SJed Brown     }
6463e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
65efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
664b9ad928SBarry Smith     while (cycles--) {
67*fcb023d4SJed Brown       ierr = PCMGMCycle_Private(pc,mglevelsin-1,transpose,reason);CHKERRQ(ierr);
684b9ad928SBarry Smith     }
6963e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
70*fcb023d4SJed Brown     if (!transpose) {
7131567311SBarry Smith       ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
72*fcb023d4SJed Brown     } else {
73*fcb023d4SJed Brown       ierr = MatInterpolateAdd(mglevels->restrct,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
74*fcb023d4SJed Brown     }
7563e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
7663e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
77*fcb023d4SJed Brown     if (!transpose) {
7831567311SBarry Smith       ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
79c0decd05SBarry Smith       ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr);
80*fcb023d4SJed Brown     } else {
81*fcb023d4SJed Brown       ierr = KSPSolveTranspose(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
82*fcb023d4SJed Brown       ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr);
83*fcb023d4SJed Brown     }
8441b6fd38SMatthew G. Knepley     if (mglevels->cr) {
8541b6fd38SMatthew G. Knepley       /* TODO Turn on copy and turn off noisy if we have an exact solution
8641b6fd38SMatthew G. Knepley       ierr = VecCopy(mglevels->x, mglevels->crx);CHKERRQ(ierr);
8741b6fd38SMatthew G. Knepley       ierr = VecCopy(mglevels->b, mglevels->crb);CHKERRQ(ierr); */
8841b6fd38SMatthew G. Knepley       ierr = KSPSetNoisy_Private(mglevels->crx);CHKERRQ(ierr);
8941b6fd38SMatthew G. Knepley       ierr = KSPSolve(mglevels->cr,mglevels->crb,mglevels->crx);CHKERRQ(ierr);    /* compatible relaxation */
9041b6fd38SMatthew G. Knepley       ierr = KSPCheckSolve(mglevels->cr,pc,mglevels->crx);CHKERRQ(ierr);
9141b6fd38SMatthew G. Knepley     }
9263e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
934b9ad928SBarry Smith   }
944b9ad928SBarry Smith   PetscFunctionReturn(0);
954b9ad928SBarry Smith }
964b9ad928SBarry Smith 
97ace3abfcSBarry 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)
984b9ad928SBarry Smith {
99f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
100f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
101dfbe8321SBarry Smith   PetscErrorCode ierr;
102391689abSStefano Zampini   PC             tpc;
103391689abSStefano Zampini   PetscBool      changeu,changed;
104f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
1054b9ad928SBarry Smith 
1064b9ad928SBarry Smith   PetscFunctionBegin;
107a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
108a762d673SBarry Smith   for (i=0; i<levels; i++) {
109a762d673SBarry Smith     if (!mglevels[i]->A) {
110a762d673SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
111a762d673SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
112a762d673SBarry Smith     }
113a762d673SBarry Smith   }
114391689abSStefano Zampini 
115391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
116391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
117391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
118391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
119391689abSStefano Zampini   if (!changed && !changeu) {
120391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
121f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
122391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
123391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
124391689abSStefano Zampini       Vec *vec;
125391689abSStefano Zampini 
126391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
127391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
1287ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
129391689abSStefano Zampini     }
130391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
131391689abSStefano Zampini   }
132f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
1334b9ad928SBarry Smith 
13431567311SBarry Smith   mg->rtol   = rtol;
13531567311SBarry Smith   mg->abstol = abstol;
13631567311SBarry Smith   mg->dtol   = dtol;
1374b9ad928SBarry Smith   if (rtol) {
1384b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
1394b9ad928SBarry Smith     PetscReal rnorm;
1407319c654SBarry Smith     if (zeroguess) {
1417319c654SBarry Smith       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
1427319c654SBarry Smith     } else {
143f3fbd535SBarry Smith       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
1444b9ad928SBarry Smith       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
1457319c654SBarry Smith     }
14631567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
1472fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
1482fa5cd67SKarl Rupp   else mg->ttol = 0.0;
1494b9ad928SBarry Smith 
1504d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
151336babb1SJed Brown      stop prematurely due to small residual */
1524d0a8057SBarry Smith   for (i=1; i<levels; i++) {
153f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
154f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
15523067569SBarry Smith       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
15623067569SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
157f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1584b9ad928SBarry Smith     }
1594d0a8057SBarry Smith   }
1604d0a8057SBarry Smith 
1614d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1624d0a8057SBarry Smith   for (i=0; i<its; i++) {
163*fcb023d4SJed Brown     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_FALSE,reason);CHKERRQ(ierr);
1644d0a8057SBarry Smith     if (*reason) break;
1654d0a8057SBarry Smith   }
1664d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1674d0a8057SBarry Smith   *outits = i;
168391689abSStefano Zampini   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
1694b9ad928SBarry Smith   PetscFunctionReturn(0);
1704b9ad928SBarry Smith }
1714b9ad928SBarry Smith 
1723aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1733aeaf226SBarry Smith {
1743aeaf226SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1753aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1763aeaf226SBarry Smith   PetscErrorCode ierr;
177f3b08a26SMatthew G. Knepley   PetscInt       i,c,n;
1783aeaf226SBarry Smith 
1793aeaf226SBarry Smith   PetscFunctionBegin;
1803aeaf226SBarry Smith   if (mglevels) {
1813aeaf226SBarry Smith     n = mglevels[0]->levels;
1823aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1833aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1843aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1853aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
18641b6fd38SMatthew G. Knepley       ierr = VecDestroy(&mglevels[i]->crx);CHKERRQ(ierr);
18741b6fd38SMatthew G. Knepley       ierr = VecDestroy(&mglevels[i]->crb);CHKERRQ(ierr);
1883aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1893aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
1900de8493bSLawrence Mitchell       ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr);
19173250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1923aeaf226SBarry Smith     }
19341b6fd38SMatthew G. Knepley     ierr = VecDestroy(&mglevels[n-1]->crx);CHKERRQ(ierr);
19441b6fd38SMatthew G. Knepley     ierr = VecDestroy(&mglevels[n-1]->crb);CHKERRQ(ierr);
195391689abSStefano Zampini     /* this is not null only if the smoother on the finest level
196391689abSStefano Zampini        changes the rhs during PreSolve */
197391689abSStefano Zampini     ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr);
1983aeaf226SBarry Smith 
1993aeaf226SBarry Smith     for (i=0; i<n; i++) {
200f3b08a26SMatthew G. Knepley       if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {ierr = VecDestroy(&mglevels[i]->coarseSpace[c]);CHKERRQ(ierr);}
201f3b08a26SMatthew G. Knepley       ierr = PetscFree(mglevels[i]->coarseSpace);CHKERRQ(ierr);
202f3b08a26SMatthew G. Knepley       mglevels[i]->coarseSpace = NULL;
2033aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
2043aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2053aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
2063aeaf226SBarry Smith       }
2073aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
20841b6fd38SMatthew G. Knepley       if (mglevels[i]->cr) {ierr = KSPReset(mglevels[i]->cr);CHKERRQ(ierr);}
2093aeaf226SBarry Smith     }
210f3b08a26SMatthew G. Knepley     mg->Nc = 0;
2113aeaf226SBarry Smith   }
2123aeaf226SBarry Smith   PetscFunctionReturn(0);
2133aeaf226SBarry Smith }
2143aeaf226SBarry Smith 
21541b6fd38SMatthew G. Knepley /* Implementing CR
21641b6fd38SMatthew G. Knepley 
21741b6fd38SMatthew 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
21841b6fd38SMatthew G. Knepley 
21941b6fd38SMatthew G. Knepley   Inj^T (Inj Inj^T)^{-1} Inj
22041b6fd38SMatthew G. Knepley 
22141b6fd38SMatthew G. Knepley and if Inj is a VecScatter, as it is now in PETSc, we have
22241b6fd38SMatthew G. Knepley 
22341b6fd38SMatthew G. Knepley   Inj^T Inj
22441b6fd38SMatthew G. Knepley 
22541b6fd38SMatthew G. Knepley and
22641b6fd38SMatthew G. Knepley 
22741b6fd38SMatthew G. Knepley   S = I - Inj^T Inj
22841b6fd38SMatthew G. Knepley 
22941b6fd38SMatthew G. Knepley since
23041b6fd38SMatthew G. Knepley 
23141b6fd38SMatthew G. Knepley   Inj S = Inj - (Inj Inj^T) Inj = 0.
23241b6fd38SMatthew G. Knepley 
23341b6fd38SMatthew G. Knepley Brannick suggests
23441b6fd38SMatthew G. Knepley 
23541b6fd38SMatthew G. Knepley   A \to S^T A S  \qquad\mathrm{and}\qquad M \to S^T M S
23641b6fd38SMatthew G. Knepley 
23741b6fd38SMatthew 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
23841b6fd38SMatthew G. Knepley 
23941b6fd38SMatthew G. Knepley   M^{-1} A \to S M^{-1} A S
24041b6fd38SMatthew G. Knepley 
24141b6fd38SMatthew G. Knepley In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left.
24241b6fd38SMatthew G. Knepley 
24341b6fd38SMatthew G. Knepley   Check: || Inj P - I ||_F < tol
24441b6fd38SMatthew G. Knepley   Check: In general, Inj Inj^T = I
24541b6fd38SMatthew G. Knepley */
24641b6fd38SMatthew G. Knepley 
24741b6fd38SMatthew G. Knepley typedef struct {
24841b6fd38SMatthew G. Knepley   PC       mg;  /* The PCMG object */
24941b6fd38SMatthew G. Knepley   PetscInt l;   /* The multigrid level for this solver */
25041b6fd38SMatthew G. Knepley   Mat      Inj; /* The injection matrix */
25141b6fd38SMatthew G. Knepley   Mat      S;   /* I - Inj^T Inj */
25241b6fd38SMatthew G. Knepley } CRContext;
25341b6fd38SMatthew G. Knepley 
25441b6fd38SMatthew G. Knepley static PetscErrorCode CRSetup_Private(PC pc)
25541b6fd38SMatthew G. Knepley {
25641b6fd38SMatthew G. Knepley   CRContext     *ctx;
25741b6fd38SMatthew G. Knepley   Mat            It;
25841b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
25941b6fd38SMatthew G. Knepley 
26041b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
26141b6fd38SMatthew G. Knepley   ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr);
26241b6fd38SMatthew G. Knepley   ierr = PCMGGetInjection(ctx->mg, ctx->l, &It);CHKERRQ(ierr);
26341b6fd38SMatthew G. Knepley   if (!It) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG");
26441b6fd38SMatthew G. Knepley   ierr = MatCreateTranspose(It, &ctx->Inj);CHKERRQ(ierr);
26541b6fd38SMatthew G. Knepley   ierr = MatCreateNormal(ctx->Inj, &ctx->S);CHKERRQ(ierr);
26641b6fd38SMatthew G. Knepley   ierr = MatScale(ctx->S, -1.0);CHKERRQ(ierr);
26741b6fd38SMatthew G. Knepley   ierr = MatShift(ctx->S,  1.0);CHKERRQ(ierr);
26841b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
26941b6fd38SMatthew G. Knepley }
27041b6fd38SMatthew G. Knepley 
27141b6fd38SMatthew G. Knepley static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y)
27241b6fd38SMatthew G. Knepley {
27341b6fd38SMatthew G. Knepley   CRContext     *ctx;
27441b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
27541b6fd38SMatthew G. Knepley 
27641b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
27741b6fd38SMatthew G. Knepley   ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr);
27841b6fd38SMatthew G. Knepley   ierr = MatMult(ctx->S, x, y);CHKERRQ(ierr);
27941b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
28041b6fd38SMatthew G. Knepley }
28141b6fd38SMatthew G. Knepley 
28241b6fd38SMatthew G. Knepley static PetscErrorCode CRDestroy_Private(PC pc)
28341b6fd38SMatthew G. Knepley {
28441b6fd38SMatthew G. Knepley   CRContext     *ctx;
28541b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
28641b6fd38SMatthew G. Knepley 
28741b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
28841b6fd38SMatthew G. Knepley   ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr);
28941b6fd38SMatthew G. Knepley   ierr = MatDestroy(&ctx->Inj);CHKERRQ(ierr);
29041b6fd38SMatthew G. Knepley   ierr = MatDestroy(&ctx->S);CHKERRQ(ierr);
29141b6fd38SMatthew G. Knepley   ierr = PetscFree(ctx);CHKERRQ(ierr);
29241b6fd38SMatthew G. Knepley   ierr = PCShellSetContext(pc, NULL);CHKERRQ(ierr);
29341b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
29441b6fd38SMatthew G. Knepley }
29541b6fd38SMatthew G. Knepley 
29641b6fd38SMatthew G. Knepley static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr)
29741b6fd38SMatthew G. Knepley {
29841b6fd38SMatthew G. Knepley   CRContext     *ctx;
29941b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
30041b6fd38SMatthew G. Knepley 
30141b6fd38SMatthew G. Knepley   PetscFunctionBeginUser;
30241b6fd38SMatthew G. Knepley   ierr = PCCreate(PetscObjectComm((PetscObject) pc), cr);CHKERRQ(ierr);
30341b6fd38SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *cr, "S (complementary projector to injection)");CHKERRQ(ierr);
30441b6fd38SMatthew G. Knepley   ierr = PetscCalloc1(1, &ctx);CHKERRQ(ierr);
30541b6fd38SMatthew G. Knepley   ctx->mg = pc;
30641b6fd38SMatthew G. Knepley   ctx->l  = l;
30741b6fd38SMatthew G. Knepley   ierr = PCSetType(*cr, PCSHELL);CHKERRQ(ierr);
30841b6fd38SMatthew G. Knepley   ierr = PCShellSetContext(*cr, ctx);CHKERRQ(ierr);
30941b6fd38SMatthew G. Knepley   ierr = PCShellSetApply(*cr, CRApply_Private);CHKERRQ(ierr);
31041b6fd38SMatthew G. Knepley   ierr = PCShellSetSetUp(*cr, CRSetup_Private);CHKERRQ(ierr);
31141b6fd38SMatthew G. Knepley   ierr = PCShellSetDestroy(*cr, CRDestroy_Private);CHKERRQ(ierr);
31241b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
31341b6fd38SMatthew G. Knepley }
31441b6fd38SMatthew G. Knepley 
31597d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
3164b9ad928SBarry Smith {
317dfbe8321SBarry Smith   PetscErrorCode ierr;
318f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
319ce94432eSBarry Smith   MPI_Comm       comm;
3203aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
32110eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
32210167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
323f3fbd535SBarry Smith   PetscInt       i;
324f3fbd535SBarry Smith   PetscMPIInt    size;
325f3fbd535SBarry Smith   const char     *prefix;
326f3fbd535SBarry Smith   PC             ipc;
3273aeaf226SBarry Smith   PetscInt       n;
3284b9ad928SBarry Smith 
3294b9ad928SBarry Smith   PetscFunctionBegin;
3300700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
331548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
332548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
3331a97d7b6SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
3343aeaf226SBarry Smith   if (mglevels) {
33510eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
3363aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
3373aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
3383aeaf226SBarry Smith     n    = mglevels[0]->levels;
3393aeaf226SBarry Smith     for (i=0; i<n; i++) {
3403aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
3413aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
3423aeaf226SBarry Smith       }
3433aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
34441b6fd38SMatthew G. Knepley       ierr = KSPDestroy(&mglevels[i]->cr);CHKERRQ(ierr);
3453aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
3463aeaf226SBarry Smith     }
3473aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
3483aeaf226SBarry Smith   }
349f3fbd535SBarry Smith 
350f3fbd535SBarry Smith   mg->nlevels = levels;
351f3fbd535SBarry Smith 
352785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
3533bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
354f3fbd535SBarry Smith 
355f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
356f3fbd535SBarry Smith 
357a9db87a2SMatthew G Knepley   mg->stageApply = 0;
358f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
359b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
3602fa5cd67SKarl Rupp 
361f3fbd535SBarry Smith     mglevels[i]->level               = i;
362f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
36310eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
364336babb1SJed Brown     mg->default_smoothu              = 2;
365336babb1SJed Brown     mg->default_smoothd              = 2;
36663e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
36763e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
36863e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
36963e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
370f3fbd535SBarry Smith 
371f3fbd535SBarry Smith     if (comms) comm = comms[i];
372d5a8f7e6SBarry Smith     if (comm != MPI_COMM_NULL) {
373f3fbd535SBarry Smith       ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
374422a814eSBarry Smith       ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
3750ee9a668SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
3760ee9a668SBarry Smith       ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
3770ee9a668SBarry Smith       ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
3780ee9a668SBarry Smith       if (i || levels == 1) {
3790ee9a668SBarry Smith         char tprefix[128];
3800ee9a668SBarry Smith 
381336babb1SJed Brown         ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
3820059c7bdSJed Brown         ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
383336babb1SJed Brown         ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
384336babb1SJed Brown         ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
385336babb1SJed Brown         ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
3860ee9a668SBarry Smith         ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
387f3fbd535SBarry Smith 
38841b6fd38SMatthew G. Knepley         ierr = PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i);CHKERRQ(ierr);
3890ee9a668SBarry Smith         ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
3900ee9a668SBarry Smith       } else {
391f3fbd535SBarry Smith         ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
392f3fbd535SBarry Smith 
3939dbfc187SHong Zhang         /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
394f3fbd535SBarry Smith         ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
395f3fbd535SBarry Smith         ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
396ffc4695bSBarry Smith         ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
397f3fbd535SBarry Smith         if (size > 1) {
398f3fbd535SBarry Smith           ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
399f3fbd535SBarry Smith         } else {
400f3fbd535SBarry Smith           ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
401f3fbd535SBarry Smith         }
402753b7fb9SBarry Smith         ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
403f3fbd535SBarry Smith       }
4043bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
405d5a8f7e6SBarry Smith     }
406f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
40731567311SBarry Smith     mg->rtol             = 0.0;
40831567311SBarry Smith     mg->abstol           = 0.0;
40931567311SBarry Smith     mg->dtol             = 0.0;
41031567311SBarry Smith     mg->ttol             = 0.0;
41131567311SBarry Smith     mg->cyclesperpcapply = 1;
412f3fbd535SBarry Smith   }
413f3fbd535SBarry Smith   mg->levels = mglevels;
41410eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
4154b9ad928SBarry Smith   PetscFunctionReturn(0);
4164b9ad928SBarry Smith }
4174b9ad928SBarry Smith 
41897d33e41SMatthew G. Knepley /*@C
41997d33e41SMatthew G. Knepley    PCMGSetLevels - Sets the number of levels to use with MG.
42097d33e41SMatthew G. Knepley    Must be called before any other MG routine.
42197d33e41SMatthew G. Knepley 
42297d33e41SMatthew G. Knepley    Logically Collective on PC
42397d33e41SMatthew G. Knepley 
42497d33e41SMatthew G. Knepley    Input Parameters:
42597d33e41SMatthew G. Knepley +  pc - the preconditioner context
42697d33e41SMatthew G. Knepley .  levels - the number of levels
42797d33e41SMatthew G. Knepley -  comms - optional communicators for each level; this is to allow solving the coarser problems
428d5a8f7e6SBarry Smith            on smaller sets of processes. For processes that are not included in the computation
429d5a8f7e6SBarry Smith            you must pass MPI_COMM_NULL.
43097d33e41SMatthew G. Knepley 
43197d33e41SMatthew G. Knepley    Level: intermediate
43297d33e41SMatthew G. Knepley 
43397d33e41SMatthew G. Knepley    Notes:
43497d33e41SMatthew G. Knepley      If the number of levels is one then the multigrid uses the -mg_levels prefix
43597d33e41SMatthew G. Knepley      for setting the level options rather than the -mg_coarse prefix.
43697d33e41SMatthew G. Knepley 
437d5a8f7e6SBarry Smith      You can free the information in comms after this routine is called.
438d5a8f7e6SBarry Smith 
439d5a8f7e6SBarry Smith      The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level
440d5a8f7e6SBarry Smith      are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on
441d5a8f7e6SBarry Smith      the two levels, rank 0 in the original communicator will pass in an array of 2 communicators
442d5a8f7e6SBarry Smith      of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators
443d5a8f7e6SBarry Smith      the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate
444d5a8f7e6SBarry Smith      in the coarse grid solve.
445d5a8f7e6SBarry Smith 
446d5a8f7e6SBarry Smith      Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one
447d5a8f7e6SBarry Smith      must take special care in providing the restriction and interpolation operation. We recommend
448d5a8f7e6SBarry Smith      providing these as two step operations; first perform a standard restriction or interpolation on
449d5a8f7e6SBarry Smith      the full number of ranks for that level and then use an MPI call to copy the resulting vector
450d5a8f7e6SBarry Smith      array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, not in both
451d5a8f7e6SBarry Smith      cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and
452d5a8f7e6SBarry Smith      recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries.
453d5a8f7e6SBarry Smith 
454d5a8f7e6SBarry Smith 
45597d33e41SMatthew G. Knepley .seealso: PCMGSetType(), PCMGGetLevels()
45697d33e41SMatthew G. Knepley @*/
45797d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
45897d33e41SMatthew G. Knepley {
45997d33e41SMatthew G. Knepley   PetscErrorCode ierr;
46097d33e41SMatthew G. Knepley 
46197d33e41SMatthew G. Knepley   PetscFunctionBegin;
46297d33e41SMatthew G. Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
46397d33e41SMatthew G. Knepley   if (comms) PetscValidPointer(comms,3);
46437b5128cSMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr);
46597d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
46697d33e41SMatthew G. Knepley }
46797d33e41SMatthew G. Knepley 
468c07bf074SBarry Smith 
469c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
470c07bf074SBarry Smith {
471c07bf074SBarry Smith   PetscErrorCode ierr;
472a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
473a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
474a06653b4SBarry Smith   PetscInt       i,n;
475c07bf074SBarry Smith 
476c07bf074SBarry Smith   PetscFunctionBegin;
477a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
478a06653b4SBarry Smith   if (mglevels) {
479a06653b4SBarry Smith     n = mglevels[0]->levels;
480a06653b4SBarry Smith     for (i=0; i<n; i++) {
481a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
4826bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
483a06653b4SBarry Smith       }
4846bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
48541b6fd38SMatthew G. Knepley       ierr = KSPDestroy(&mglevels[i]->cr);CHKERRQ(ierr);
486a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
487a06653b4SBarry Smith     }
488a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
489a06653b4SBarry Smith   }
490c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
491fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);CHKERRQ(ierr);
492fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);CHKERRQ(ierr);
493f3fbd535SBarry Smith   PetscFunctionReturn(0);
494f3fbd535SBarry Smith }
495f3fbd535SBarry Smith 
496f3fbd535SBarry Smith 
497f3fbd535SBarry Smith /*
498f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
499f3fbd535SBarry Smith              or full cycle of multigrid.
500f3fbd535SBarry Smith 
501f3fbd535SBarry Smith   Note:
502f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
503f3fbd535SBarry Smith */
504*fcb023d4SJed Brown static PetscErrorCode PCApply_MG_Internal(PC pc,Vec b,Vec x,PetscBool transpose)
505f3fbd535SBarry Smith {
506f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
507f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
508f3fbd535SBarry Smith   PetscErrorCode ierr;
509391689abSStefano Zampini   PC             tpc;
510f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
511391689abSStefano Zampini   PetscBool      changeu,changed;
512f3fbd535SBarry Smith 
513f3fbd535SBarry Smith   PetscFunctionBegin;
514a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
515e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
516298cc208SBarry Smith   for (i=0; i<levels; i++) {
517e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
51823ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
519298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
520e1d8e5deSBarry Smith     }
521e1d8e5deSBarry Smith   }
522e1d8e5deSBarry Smith 
523391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
524391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
525391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
526391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
527391689abSStefano Zampini   if (!changeu && !changed) {
528391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
529f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
530391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
531391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
532391689abSStefano Zampini       Vec *vec;
533391689abSStefano Zampini 
534391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
535391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
5367ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
537391689abSStefano Zampini     }
538391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
539391689abSStefano Zampini   }
540f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
541391689abSStefano Zampini 
54231567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
543*fcb023d4SJed Brown     ierr = VecZeroEntries(x);CHKERRQ(ierr);
54431567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
545*fcb023d4SJed Brown       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,transpose,NULL);CHKERRQ(ierr);
546f3fbd535SBarry Smith     }
5472fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
548*fcb023d4SJed Brown     ierr = PCMGACycle_Private(pc,mglevels,transpose);CHKERRQ(ierr);
5492fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
550*fcb023d4SJed Brown     ierr = PCMGKCycle_Private(pc,mglevels,transpose);CHKERRQ(ierr);
5512fa5cd67SKarl Rupp   } else {
552*fcb023d4SJed Brown     ierr = PCMGFCycle_Private(pc,mglevels,transpose);CHKERRQ(ierr);
553f3fbd535SBarry Smith   }
554a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
555391689abSStefano Zampini   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
556f3fbd535SBarry Smith   PetscFunctionReturn(0);
557f3fbd535SBarry Smith }
558f3fbd535SBarry Smith 
559*fcb023d4SJed Brown static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
560*fcb023d4SJed Brown {
561*fcb023d4SJed Brown     PetscErrorCode ierr;
562*fcb023d4SJed Brown 
563*fcb023d4SJed Brown     PetscFunctionBegin;
564*fcb023d4SJed Brown     ierr = PCApply_MG_Internal(pc,b,x,PETSC_FALSE);CHKERRQ(ierr);
565*fcb023d4SJed Brown     PetscFunctionReturn(0);
566*fcb023d4SJed Brown }
567*fcb023d4SJed Brown 
568*fcb023d4SJed Brown static PetscErrorCode PCApplyTranspose_MG(PC pc,Vec b,Vec x)
569*fcb023d4SJed Brown {
570*fcb023d4SJed Brown     PetscErrorCode ierr;
571*fcb023d4SJed Brown 
572*fcb023d4SJed Brown     PetscFunctionBegin;
573*fcb023d4SJed Brown     ierr = PCApply_MG_Internal(pc,b,x,PETSC_TRUE);CHKERRQ(ierr);
574*fcb023d4SJed Brown     PetscFunctionReturn(0);
575*fcb023d4SJed Brown }
576f3fbd535SBarry Smith 
5774416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
578f3fbd535SBarry Smith {
579f3fbd535SBarry Smith   PetscErrorCode   ierr;
580710315b6SLawrence Mitchell   PetscInt         levels,cycles;
581f3b08a26SMatthew G. Knepley   PetscBool        flg, flg2;
582f3fbd535SBarry Smith   PC_MG            *mg = (PC_MG*)pc->data;
5833d3eaba7SBarry Smith   PC_MG_Levels     **mglevels;
584f3fbd535SBarry Smith   PCMGType         mgtype;
585f3fbd535SBarry Smith   PCMGCycleType    mgctype;
5862134b1e4SBarry Smith   PCMGGalerkinType gtype;
587f3fbd535SBarry Smith 
588f3fbd535SBarry Smith   PetscFunctionBegin;
5891d743356SStefano Zampini   levels = PetscMax(mg->nlevels,1);
590e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
591f3fbd535SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
5921a97d7b6SStefano Zampini   if (!flg && !mg->levels && pc->dm) {
593eb3f98d2SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
594eb3f98d2SBarry Smith     levels++;
5953aeaf226SBarry Smith     mg->usedmfornumberoflevels = PETSC_TRUE;
596eb3f98d2SBarry Smith   }
5970298fd71SBarry Smith   ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5983aeaf226SBarry Smith   mglevels = mg->levels;
5993aeaf226SBarry Smith 
600f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
601f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
602f3fbd535SBarry Smith   if (flg) {
603f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
6042fa5cd67SKarl Rupp   }
6052134b1e4SBarry Smith   gtype = mg->galerkin;
6062134b1e4SBarry Smith   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
6072134b1e4SBarry Smith   if (flg) {
6082134b1e4SBarry Smith     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
609f3fbd535SBarry Smith   }
610f3b08a26SMatthew G. Knepley   flg2 = PETSC_FALSE;
611f3b08a26SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr);
612f3b08a26SMatthew G. Knepley   if (flg) {ierr = PCMGSetAdaptInterpolation(pc, flg2);CHKERRQ(ierr);}
613f3b08a26SMatthew 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);
614f3b08a26SMatthew 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);
615f3b08a26SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);CHKERRQ(ierr);
61641b6fd38SMatthew G. Knepley   flg2 = PETSC_FALSE;
61741b6fd38SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr);
61841b6fd38SMatthew G. Knepley   if (flg) {ierr = PCMGSetAdaptCR(pc, flg2);CHKERRQ(ierr);}
619f56b1265SBarry Smith   flg = PETSC_FALSE;
6208e5aa403SBarry Smith   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
621f442ab6aSBarry Smith   if (flg) {
622f442ab6aSBarry Smith     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
623f442ab6aSBarry Smith   }
62431567311SBarry Smith   mgtype = mg->am;
625f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
626f3fbd535SBarry Smith   if (flg) {
627f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
628f3fbd535SBarry Smith   }
62931567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
6305363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
631f3fbd535SBarry Smith     if (flg) {
632f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
633f3fbd535SBarry Smith     }
634f3fbd535SBarry Smith   }
635f3fbd535SBarry Smith   flg  = PETSC_FALSE;
6360298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
637f3fbd535SBarry Smith   if (flg) {
638f3fbd535SBarry Smith     PetscInt i;
639f3fbd535SBarry Smith     char     eventname[128];
6401a97d7b6SStefano Zampini 
641f3fbd535SBarry Smith     levels = mglevels[0]->levels;
642f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
643f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
64463e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
645f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
64663e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
647f3fbd535SBarry Smith       if (i) {
648f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
64963e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
650f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
65163e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
652f3fbd535SBarry Smith       }
653f3fbd535SBarry Smith     }
654eec35531SMatthew G Knepley 
655ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
656eec35531SMatthew G Knepley     {
657eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
658eec35531SMatthew G Knepley       PetscStageLog stageLog;
659eec35531SMatthew G Knepley       PetscInt      st;
660eec35531SMatthew G Knepley 
661eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
662eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
663eec35531SMatthew G Knepley         PetscBool same;
664eec35531SMatthew G Knepley 
665eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
6662fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
667eec35531SMatthew G Knepley       }
668eec35531SMatthew G Knepley       if (!mg->stageApply) {
669eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
670eec35531SMatthew G Knepley       }
671eec35531SMatthew G Knepley     }
672ec60ca73SSatish Balay #endif
673f3fbd535SBarry Smith   }
674f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
675f3b08a26SMatthew G. Knepley   /* Check option consistency */
676f3b08a26SMatthew G. Knepley   ierr = PCMGGetGalerkin(pc, &gtype);CHKERRQ(ierr);
677f3b08a26SMatthew G. Knepley   ierr = PCMGGetAdaptInterpolation(pc, &flg);CHKERRQ(ierr);
678f3b08a26SMatthew 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");
679f3fbd535SBarry Smith   PetscFunctionReturn(0);
680f3fbd535SBarry Smith }
681f3fbd535SBarry Smith 
6820a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL};
6830a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL};
6840a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL};
685f3b08a26SMatthew G. Knepley const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL};
686f3fbd535SBarry Smith 
6879804daf3SBarry Smith #include <petscdraw.h>
688f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
689f3fbd535SBarry Smith {
690f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
691f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
692f3fbd535SBarry Smith   PetscErrorCode ierr;
693e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
694219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
695f3fbd535SBarry Smith 
696f3fbd535SBarry Smith   PetscFunctionBegin;
697251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
6985b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
699219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
700f3fbd535SBarry Smith   if (iascii) {
701e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
702efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
70331567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
70431567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
705f3fbd535SBarry Smith     }
7062134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
707f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
7082134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
7092134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
7102134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
7112134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
7122134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
7132134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
7144f66f45eSBarry Smith     } else {
7154f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
716f3fbd535SBarry Smith     }
7175adeb434SBarry Smith     if (mg->view){
7185adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
7195adeb434SBarry Smith     }
720f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
721f3fbd535SBarry Smith       if (!i) {
72289cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
723f3fbd535SBarry Smith       } else {
72489cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
725f3fbd535SBarry Smith       }
726f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
727f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
728f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
729f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
730f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
731f3fbd535SBarry Smith       } else if (i) {
73289cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
733f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
734f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
735f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
736f3fbd535SBarry Smith       }
73741b6fd38SMatthew G. Knepley       if (i && mglevels[i]->cr) {
73841b6fd38SMatthew G. Knepley         ierr = PetscViewerASCIIPrintf(viewer,"CR solver on level %D -------------------------------\n",i);CHKERRQ(ierr);
73941b6fd38SMatthew G. Knepley         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
74041b6fd38SMatthew G. Knepley         ierr = KSPView(mglevels[i]->cr,viewer);CHKERRQ(ierr);
74141b6fd38SMatthew G. Knepley         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
74241b6fd38SMatthew G. Knepley       }
743f3fbd535SBarry Smith     }
7445b0b0462SBarry Smith   } else if (isbinary) {
7455b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
7465b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
7475b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
7485b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
7495b0b0462SBarry Smith       }
7505b0b0462SBarry Smith     }
751219b1045SBarry Smith   } else if (isdraw) {
752219b1045SBarry Smith     PetscDraw draw;
753b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
754219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
755219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
7560298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
757219b1045SBarry Smith     bottom = y - th;
758219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
759b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
760219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
761219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
762219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
763b4375e8dSPeter Brune       } else {
764b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
765b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
766b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
767b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
768b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
769b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
770b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
771b4375e8dSPeter Brune       }
7720298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
7731cd381d6SBarry Smith       bottom -= th;
774219b1045SBarry Smith     }
7755b0b0462SBarry Smith   }
776f3fbd535SBarry Smith   PetscFunctionReturn(0);
777f3fbd535SBarry Smith }
778f3fbd535SBarry Smith 
779af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
780cab2e9ccSBarry Smith 
781f3fbd535SBarry Smith /*
782f3fbd535SBarry Smith     Calls setup for the KSP on each level
783f3fbd535SBarry Smith */
784f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
785f3fbd535SBarry Smith {
786f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
787f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
788f3fbd535SBarry Smith   PetscErrorCode ierr;
7891a97d7b6SStefano Zampini   PetscInt       i,n;
79098e04984SBarry Smith   PC             cpc;
7913db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
792f3fbd535SBarry Smith   Mat            dA,dB;
793f3fbd535SBarry Smith   Vec            tvec;
794218a07d4SBarry Smith   DM             *dms;
7950a545947SLisandro Dalcin   PetscViewer    viewer = NULL;
79641b6fd38SMatthew G. Knepley   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE;
797f3fbd535SBarry Smith 
798f3fbd535SBarry Smith   PetscFunctionBegin;
7991a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
8001a97d7b6SStefano Zampini   n = mglevels[0]->levels;
80101bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
8023aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
8033aeaf226SBarry Smith     PetscInt levels;
8043aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
8053aeaf226SBarry Smith     levels++;
8063aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
8070298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
8083aeaf226SBarry Smith       n        = levels;
8093aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
8103aeaf226SBarry Smith       mglevels = mg->levels;
8113aeaf226SBarry Smith     }
8123aeaf226SBarry Smith   }
81398e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
8143aeaf226SBarry Smith 
815f3fbd535SBarry Smith 
816f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
817f3fbd535SBarry Smith   /* so use those from global PC */
818f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
8190298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
82098e04984SBarry Smith   if (opsset) {
82198e04984SBarry Smith     Mat mmat;
82223ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
82398e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
82498e04984SBarry Smith   }
8254a5f07a5SMark F. Adams 
82641b6fd38SMatthew G. Knepley   /* Create CR solvers */
82741b6fd38SMatthew G. Knepley   ierr = PCMGGetAdaptCR(pc, &doCR);CHKERRQ(ierr);
82841b6fd38SMatthew G. Knepley   if (doCR) {
82941b6fd38SMatthew G. Knepley     const char *prefix;
83041b6fd38SMatthew G. Knepley 
83141b6fd38SMatthew G. Knepley     ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr);
83241b6fd38SMatthew G. Knepley     for (i = 1; i < n; ++i) {
83341b6fd38SMatthew G. Knepley       PC   ipc, cr;
83441b6fd38SMatthew G. Knepley       char crprefix[128];
83541b6fd38SMatthew G. Knepley 
83641b6fd38SMatthew G. Knepley       ierr = KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr);CHKERRQ(ierr);
83741b6fd38SMatthew G. Knepley       ierr = KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE);CHKERRQ(ierr);
83841b6fd38SMatthew G. Knepley       ierr = PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i);CHKERRQ(ierr);
83941b6fd38SMatthew G. Knepley       ierr = KSPSetOptionsPrefix(mglevels[i]->cr, prefix);CHKERRQ(ierr);
84041b6fd38SMatthew G. Knepley       ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
84141b6fd38SMatthew G. Knepley       ierr = KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV);CHKERRQ(ierr);
84241b6fd38SMatthew G. Knepley       ierr = KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL);CHKERRQ(ierr);
84341b6fd38SMatthew G. Knepley       ierr = KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED);CHKERRQ(ierr);
84441b6fd38SMatthew G. Knepley       ierr = KSPGetPC(mglevels[i]->cr, &ipc);CHKERRQ(ierr);
84541b6fd38SMatthew G. Knepley 
84641b6fd38SMatthew G. Knepley       ierr = PCSetType(ipc, PCCOMPOSITE);CHKERRQ(ierr);
84741b6fd38SMatthew G. Knepley       ierr = PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr);
84841b6fd38SMatthew G. Knepley       ierr = PCCompositeAddPCType(ipc, PCSOR);CHKERRQ(ierr);
84941b6fd38SMatthew G. Knepley       ierr = CreateCR_Private(pc, i, &cr);CHKERRQ(ierr);
85041b6fd38SMatthew G. Knepley       ierr = PCCompositeAddPC(ipc, cr);CHKERRQ(ierr);
85141b6fd38SMatthew G. Knepley       ierr = PCDestroy(&cr);CHKERRQ(ierr);
85241b6fd38SMatthew G. Knepley 
85341b6fd38SMatthew G. Knepley       ierr = KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
85441b6fd38SMatthew G. Knepley       ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE);CHKERRQ(ierr);
85541b6fd38SMatthew G. Knepley       ierr = PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i);CHKERRQ(ierr);
85641b6fd38SMatthew G. Knepley       ierr = KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix);CHKERRQ(ierr);
85741b6fd38SMatthew G. Knepley       ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr);CHKERRQ(ierr);
85841b6fd38SMatthew G. Knepley     }
85941b6fd38SMatthew G. Knepley   }
86041b6fd38SMatthew G. Knepley 
8614a5f07a5SMark F. Adams   if (!opsset) {
86271b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
8634a5f07a5SMark F. Adams     if (use_amat) {
864f3fbd535SBarry 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);
86523ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
86669858f1bSStefano Zampini     } else {
8674a5f07a5SMark 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);
86823ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
8694a5f07a5SMark F. Adams     }
8704a5f07a5SMark F. Adams   }
871f3fbd535SBarry Smith 
8729c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
8739c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
8749c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
8759c7ffb39SBarry Smith       continue;
8769c7ffb39SBarry Smith     }
8779c7ffb39SBarry Smith   }
8782134b1e4SBarry Smith 
8792134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
8802134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
8812134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
8822134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
8832134b1e4SBarry Smith   }
8842134b1e4SBarry Smith 
8852134b1e4SBarry Smith 
8869c7ffb39SBarry Smith   /*
8879c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
8889c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
8899c7ffb39SBarry Smith   */
8902134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
8912d2b81a6SBarry Smith         /* construct the interpolation from the DMs */
8922e499ae9SBarry Smith     Mat p;
89373250ac0SBarry Smith     Vec rscale;
894785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
895218a07d4SBarry Smith     dms[n-1] = pc->dm;
896ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
897ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
89841fce8fdSHong Zhang         /*
89941fce8fdSHong Zhang            Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
90041fce8fdSHong Zhang            Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
90141fce8fdSHong Zhang            But it is safe to use -dm_mat_type.
90241fce8fdSHong Zhang 
90341fce8fdSHong Zhang            The mat type should not be hardcoded like this, we need to find a better way.
90441fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
90541fce8fdSHong Zhang     */
906218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
907942e3340SBarry Smith       DMKSP     kdm;
908eab52d2bSLawrence Mitchell       PetscBool dmhasrestrict, dmhasinject;
9093c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
9102134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
911c27ee7a3SPatrick Farrell       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
912c27ee7a3SPatrick Farrell         ierr = KSPSetDM(mglevels[i]->smoothu,dms[i]);CHKERRQ(ierr);
913c27ee7a3SPatrick Farrell         if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);CHKERRQ(ierr);}
914c27ee7a3SPatrick Farrell       }
91541b6fd38SMatthew G. Knepley       if (mglevels[i]->cr) {
91641b6fd38SMatthew G. Knepley         ierr = KSPSetDM(mglevels[i]->cr,dms[i]);CHKERRQ(ierr);
91741b6fd38SMatthew G. Knepley         if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE);CHKERRQ(ierr);}
91841b6fd38SMatthew G. Knepley       }
919942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
920d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
921d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
9220298fd71SBarry Smith       kdm->ops->computerhs = NULL;
9230298fd71SBarry Smith       kdm->rhsctx          = NULL;
92424c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
925e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
9262d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
92700ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
92873250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
9296bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
930218a07d4SBarry Smith       }
9313ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
9323ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
9333ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
9343ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
9353ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
9363ad4599aSBarry Smith       }
937eab52d2bSLawrence Mitchell       ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr);
938eab52d2bSLawrence Mitchell       if (dmhasinject && !mglevels[i+1]->inject){
939eab52d2bSLawrence Mitchell         ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr);
940eab52d2bSLawrence Mitchell         ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr);
941eab52d2bSLawrence Mitchell         ierr = MatDestroy(&p);CHKERRQ(ierr);
942eab52d2bSLawrence Mitchell       }
94324c3aa18SBarry Smith     }
9442d2b81a6SBarry Smith 
945ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
9462d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
947ce4cda84SJed Brown   }
948cab2e9ccSBarry Smith 
949ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
9502134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
951cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
952cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
953c27ee7a3SPatrick Farrell     if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) {
954c27ee7a3SPatrick Farrell       ierr = KSPSetDM(mglevels[n-1]->smoothu,pc->dm);CHKERRQ(ierr);
955c27ee7a3SPatrick Farrell       ierr = KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);CHKERRQ(ierr);
956c27ee7a3SPatrick Farrell     }
95741b6fd38SMatthew G. Knepley     if (mglevels[n-1]->cr) {
95841b6fd38SMatthew G. Knepley       ierr = KSPSetDM(mglevels[n-1]->cr,pc->dm);CHKERRQ(ierr);
95941b6fd38SMatthew G. Knepley       ierr = KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE);CHKERRQ(ierr);
96041b6fd38SMatthew G. Knepley     }
961218a07d4SBarry Smith   }
962218a07d4SBarry Smith 
9632134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
9642134b1e4SBarry Smith     Mat       A,B;
9652134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
9662134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
9672134b1e4SBarry Smith 
9682134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
9692134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
9702134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
971f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
972b9367914SBarry 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");
973b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
974b9367914SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
975b9367914SBarry Smith       }
976b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
977b9367914SBarry Smith         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
978b9367914SBarry Smith       }
9792134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
9802134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
9812134b1e4SBarry Smith       }
9822134b1e4SBarry Smith       if (doA) {
9832df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
9842134b1e4SBarry Smith       }
9852134b1e4SBarry Smith       if (doB) {
9862df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
987b9367914SBarry Smith       }
9882134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
9892134b1e4SBarry Smith       if (!doA && dAeqdB) {
9902134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
9912134b1e4SBarry Smith         A = B;
9922134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
9932134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
9942134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
995b9367914SBarry Smith       }
9962134b1e4SBarry Smith       if (!doB && dAeqdB) {
9972134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
9982134b1e4SBarry Smith         B = A;
9992134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
100023ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
10012134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
10027d537d92SJed Brown       }
10032134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
10042134b1e4SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
10052134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
10062134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
10072134b1e4SBarry Smith       }
10082134b1e4SBarry Smith       dA = A;
1009f3fbd535SBarry Smith       dB = B;
1010f3fbd535SBarry Smith     }
1011f3fbd535SBarry Smith   }
1012f3b08a26SMatthew G. Knepley 
1013f3b08a26SMatthew G. Knepley 
1014f3b08a26SMatthew G. Knepley   /* Adapt interpolation matrices */
1015f3b08a26SMatthew G. Knepley   if (mg->adaptInterpolation) {
1016f3b08a26SMatthew G. Knepley     mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */
1017f3b08a26SMatthew G. Knepley     for (i = 0; i < n; ++i) {
1018f3b08a26SMatthew G. Knepley       ierr = PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);CHKERRQ(ierr);
1019f3b08a26SMatthew 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);}
1020f3b08a26SMatthew G. Knepley     }
1021f3b08a26SMatthew G. Knepley     for (i = n-2; i > -1; --i) {
1022f3b08a26SMatthew G. Knepley       ierr = PCMGRecomputeLevelOperators_Internal(pc, i);CHKERRQ(ierr);
1023f3b08a26SMatthew G. Knepley     }
1024f3b08a26SMatthew G. Knepley   }
1025f3b08a26SMatthew G. Knepley 
10262134b1e4SBarry Smith   if (needRestricts && pc->dm) {
1027caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
1028ccceb50cSJed Brown       DM  dmfine,dmcoarse;
1029caa4e7f2SJed Brown       Mat Restrict,Inject;
1030caa4e7f2SJed Brown       Vec rscale;
1031ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
1032ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
1033caa4e7f2SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
1034caa4e7f2SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
1035eab52d2bSLawrence Mitchell       ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr);
1036caa4e7f2SJed Brown       ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
1037caa4e7f2SJed Brown     }
1038caa4e7f2SJed Brown   }
1039f3fbd535SBarry Smith 
1040f3fbd535SBarry Smith   if (!pc->setupcalled) {
1041f3fbd535SBarry Smith     for (i=0; i<n; i++) {
1042f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
1043f3fbd535SBarry Smith     }
1044f3fbd535SBarry Smith     for (i=1; i<n; i++) {
1045f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
1046f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
1047f3fbd535SBarry Smith       }
104841b6fd38SMatthew G. Knepley       if (mglevels[i]->cr) {
104941b6fd38SMatthew G. Knepley         ierr = KSPSetFromOptions(mglevels[i]->cr);CHKERRQ(ierr);
105041b6fd38SMatthew G. Knepley       }
1051f3fbd535SBarry Smith     }
10523ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
1053f3fbd535SBarry Smith     for (i=1; i<n; i++) {
10543ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
10553ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
1056f3fbd535SBarry Smith     }
1057f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
1058f3fbd535SBarry Smith       if (!mglevels[i]->b) {
1059f3fbd535SBarry Smith         Vec *vec;
10602a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
1061f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
10626bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
1063f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
1064f3fbd535SBarry Smith       }
1065f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
1066f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
1067f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
10686bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
1069f3fbd535SBarry Smith       }
1070f3fbd535SBarry Smith       if (!mglevels[i]->x) {
1071f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
1072f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
10736bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
1074f3fbd535SBarry Smith       }
107541b6fd38SMatthew G. Knepley       if (doCR) {
107641b6fd38SMatthew G. Knepley         ierr = VecDuplicate(mglevels[i]->b,&mglevels[i]->crx);CHKERRQ(ierr);
107741b6fd38SMatthew G. Knepley         ierr = VecDuplicate(mglevels[i]->b,&mglevels[i]->crb);CHKERRQ(ierr);
107841b6fd38SMatthew G. Knepley         ierr = VecZeroEntries(mglevels[i]->crb);CHKERRQ(ierr);
107941b6fd38SMatthew G. Knepley       }
1080f3fbd535SBarry Smith     }
1081f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
1082f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
1083f3fbd535SBarry Smith       Vec *vec;
10842a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
1085f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
10866bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
1087f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
1088f3fbd535SBarry Smith     }
108941b6fd38SMatthew G. Knepley     if (doCR) {
109041b6fd38SMatthew G. Knepley       ierr = VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx);CHKERRQ(ierr);
109141b6fd38SMatthew G. Knepley       ierr = VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb);CHKERRQ(ierr);
109241b6fd38SMatthew G. Knepley       ierr = VecZeroEntries(mglevels[n-1]->crb);CHKERRQ(ierr);
109341b6fd38SMatthew G. Knepley     }
1094f3fbd535SBarry Smith   }
1095f3fbd535SBarry Smith 
109698e04984SBarry Smith   if (pc->dm) {
109798e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
109898e04984SBarry Smith     for (i=0; i<n-1; i++) {
109998e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
110098e04984SBarry Smith     }
110198e04984SBarry Smith   }
1102f3fbd535SBarry Smith 
1103f3fbd535SBarry Smith   for (i=1; i<n; i++) {
11042a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
1105f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
1106f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
1107f3fbd535SBarry Smith     }
110841b6fd38SMatthew G. Knepley     if (mglevels[i]->cr) {ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);CHKERRQ(ierr);}
110963e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1110f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
1111c0decd05SBarry Smith     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1112899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
1113899639b0SHong Zhang     }
111463e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1115d42688cbSBarry Smith     if (!mglevels[i]->residual) {
1116d42688cbSBarry Smith       Mat mat;
111713842ffbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
111854b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
1119d42688cbSBarry Smith     }
1120*fcb023d4SJed Brown     if (!mglevels[i]->residualtranspose) {
1121*fcb023d4SJed Brown       Mat mat;
1122*fcb023d4SJed Brown       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
1123*fcb023d4SJed Brown       ierr = PCMGSetResidualTranspose(pc,i,PCMGResidualTransposeDefault,mat);CHKERRQ(ierr);
1124*fcb023d4SJed Brown     }
1125f3fbd535SBarry Smith   }
1126f3fbd535SBarry Smith   for (i=1; i<n; i++) {
1127f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
1128f3fbd535SBarry Smith       Mat downmat,downpmat;
1129f3fbd535SBarry Smith 
1130f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
11310298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
1132f3fbd535SBarry Smith       if (!opsset) {
113323ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
113423ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
1135f3fbd535SBarry Smith       }
1136f3fbd535SBarry Smith 
1137f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
113863e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1139f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
1140c0decd05SBarry Smith       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
1141899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
1142899639b0SHong Zhang       }
114363e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1144f3fbd535SBarry Smith     }
114541b6fd38SMatthew G. Knepley     if (mglevels[i]->cr) {
114641b6fd38SMatthew G. Knepley       Mat downmat,downpmat;
114741b6fd38SMatthew G. Knepley 
114841b6fd38SMatthew G. Knepley       /* check if operators have been set for up, if not use down operators to set them */
114941b6fd38SMatthew G. Knepley       ierr = KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL);CHKERRQ(ierr);
115041b6fd38SMatthew G. Knepley       if (!opsset) {
115141b6fd38SMatthew G. Knepley         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
115241b6fd38SMatthew G. Knepley         ierr = KSPSetOperators(mglevels[i]->cr,downmat,downpmat);CHKERRQ(ierr);
115341b6fd38SMatthew G. Knepley       }
115441b6fd38SMatthew G. Knepley 
115541b6fd38SMatthew G. Knepley       ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);CHKERRQ(ierr);
115641b6fd38SMatthew G. Knepley       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
115741b6fd38SMatthew G. Knepley       ierr = KSPSetUp(mglevels[i]->cr);CHKERRQ(ierr);
115841b6fd38SMatthew G. Knepley       if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) {
115941b6fd38SMatthew G. Knepley         pc->failedreason = PC_SUBPC_ERROR;
116041b6fd38SMatthew G. Knepley       }
116141b6fd38SMatthew G. Knepley       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
116241b6fd38SMatthew G. Knepley     }
1163f3fbd535SBarry Smith   }
1164f3fbd535SBarry Smith 
116563e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1166f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
1167c0decd05SBarry Smith   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1168899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
1169899639b0SHong Zhang   }
117063e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1171f3fbd535SBarry Smith 
1172f3fbd535SBarry Smith   /*
1173f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
1174e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
1175f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
1176f3fbd535SBarry Smith 
1177f3fbd535SBarry Smith    Only support one or the other at the same time.
1178f3fbd535SBarry Smith   */
1179f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
1180c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
1181ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
1182f3fbd535SBarry Smith   dump = PETSC_FALSE;
1183f3fbd535SBarry Smith #endif
1184c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
1185ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
1186f3fbd535SBarry Smith 
1187f3fbd535SBarry Smith   if (viewer) {
1188f3fbd535SBarry Smith     for (i=1; i<n; i++) {
1189f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
1190f3fbd535SBarry Smith     }
1191f3fbd535SBarry Smith     for (i=0; i<n; i++) {
1192f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
1193f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
1194f3fbd535SBarry Smith     }
1195f3fbd535SBarry Smith   }
1196f3fbd535SBarry Smith   PetscFunctionReturn(0);
1197f3fbd535SBarry Smith }
1198f3fbd535SBarry Smith 
1199f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
1200f3fbd535SBarry Smith 
120197d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
120297d33e41SMatthew G. Knepley {
120397d33e41SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
120497d33e41SMatthew G. Knepley 
120597d33e41SMatthew G. Knepley   PetscFunctionBegin;
120697d33e41SMatthew G. Knepley   *levels = mg->nlevels;
120797d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
120897d33e41SMatthew G. Knepley }
120997d33e41SMatthew G. Knepley 
12104b9ad928SBarry Smith /*@
121197177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
12124b9ad928SBarry Smith 
12134b9ad928SBarry Smith    Not Collective
12144b9ad928SBarry Smith 
12154b9ad928SBarry Smith    Input Parameter:
12164b9ad928SBarry Smith .  pc - the preconditioner context
12174b9ad928SBarry Smith 
12184b9ad928SBarry Smith    Output parameter:
12194b9ad928SBarry Smith .  levels - the number of levels
12204b9ad928SBarry Smith 
12214b9ad928SBarry Smith    Level: advanced
12224b9ad928SBarry Smith 
122397177400SBarry Smith .seealso: PCMGSetLevels()
12244b9ad928SBarry Smith @*/
12257087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
12264b9ad928SBarry Smith {
1227e952c529SMatthew G. Knepley   PetscErrorCode ierr;
12284b9ad928SBarry Smith 
12294b9ad928SBarry Smith   PetscFunctionBegin;
12300700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
12314482741eSBarry Smith   PetscValidIntPointer(levels,2);
123297d33e41SMatthew G. Knepley   *levels = 0;
123337b5128cSMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr);
12344b9ad928SBarry Smith   PetscFunctionReturn(0);
12354b9ad928SBarry Smith }
12364b9ad928SBarry Smith 
12374b9ad928SBarry Smith /*@
123897177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
12394b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
12404b9ad928SBarry Smith 
1241ad4df100SBarry Smith    Logically Collective on PC
12424b9ad928SBarry Smith 
12434b9ad928SBarry Smith    Input Parameters:
12444b9ad928SBarry Smith +  pc - the preconditioner context
12459dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
12469dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
12474b9ad928SBarry Smith 
12484b9ad928SBarry Smith    Options Database Key:
12494b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
12504b9ad928SBarry Smith    additive, full, kaskade
12514b9ad928SBarry Smith 
12524b9ad928SBarry Smith    Level: advanced
12534b9ad928SBarry Smith 
125497177400SBarry Smith .seealso: PCMGSetLevels()
12554b9ad928SBarry Smith @*/
12567087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
12574b9ad928SBarry Smith {
1258f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
12594b9ad928SBarry Smith 
12604b9ad928SBarry Smith   PetscFunctionBegin;
12610700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1262c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
126331567311SBarry Smith   mg->am = form;
12649dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1265c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
1266c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1267c60c7ad4SBarry Smith }
1268c60c7ad4SBarry Smith 
1269c60c7ad4SBarry Smith /*@
1270c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
1271c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
1272c60c7ad4SBarry Smith 
1273c60c7ad4SBarry Smith    Logically Collective on PC
1274c60c7ad4SBarry Smith 
1275c60c7ad4SBarry Smith    Input Parameter:
1276c60c7ad4SBarry Smith .  pc - the preconditioner context
1277c60c7ad4SBarry Smith 
1278c60c7ad4SBarry Smith    Output Parameter:
1279c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
1280c60c7ad4SBarry Smith 
1281c60c7ad4SBarry Smith 
1282c60c7ad4SBarry Smith    Level: advanced
1283c60c7ad4SBarry Smith 
1284c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
1285c60c7ad4SBarry Smith @*/
1286c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1287c60c7ad4SBarry Smith {
1288c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1289c60c7ad4SBarry Smith 
1290c60c7ad4SBarry Smith   PetscFunctionBegin;
1291c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1292c60c7ad4SBarry Smith   *type = mg->am;
12934b9ad928SBarry Smith   PetscFunctionReturn(0);
12944b9ad928SBarry Smith }
12954b9ad928SBarry Smith 
12964b9ad928SBarry Smith /*@
12970d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
12984b9ad928SBarry Smith    complicated cycling.
12994b9ad928SBarry Smith 
1300ad4df100SBarry Smith    Logically Collective on PC
13014b9ad928SBarry Smith 
13024b9ad928SBarry Smith    Input Parameters:
1303c2be2410SBarry Smith +  pc - the multigrid context
1304c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
13054b9ad928SBarry Smith 
13064b9ad928SBarry Smith    Options Database Key:
1307c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
13084b9ad928SBarry Smith 
13094b9ad928SBarry Smith    Level: advanced
13104b9ad928SBarry Smith 
13110d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
13124b9ad928SBarry Smith @*/
13137087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
13144b9ad928SBarry Smith {
1315f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
1316f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
131779416396SBarry Smith   PetscInt     i,levels;
13184b9ad928SBarry Smith 
13194b9ad928SBarry Smith   PetscFunctionBegin;
13200700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
132169fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
13221a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1323f3fbd535SBarry Smith   levels = mglevels[0]->levels;
13242fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
13254b9ad928SBarry Smith   PetscFunctionReturn(0);
13264b9ad928SBarry Smith }
13274b9ad928SBarry Smith 
13288cc2d5dfSBarry Smith /*@
13298cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
13308cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
13318cc2d5dfSBarry Smith 
1332ad4df100SBarry Smith    Logically Collective on PC
13338cc2d5dfSBarry Smith 
13348cc2d5dfSBarry Smith    Input Parameters:
13358cc2d5dfSBarry Smith +  pc - the multigrid context
13368cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
13378cc2d5dfSBarry Smith 
13388cc2d5dfSBarry Smith    Options Database Key:
1339e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
13408cc2d5dfSBarry Smith 
13418cc2d5dfSBarry Smith    Level: advanced
13428cc2d5dfSBarry Smith 
134395452b02SPatrick Sanan    Notes:
134495452b02SPatrick Sanan     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
13458cc2d5dfSBarry Smith 
13468cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
13478cc2d5dfSBarry Smith @*/
13487087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
13498cc2d5dfSBarry Smith {
1350f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
13518cc2d5dfSBarry Smith 
13528cc2d5dfSBarry Smith   PetscFunctionBegin;
13530700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1354c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
13552a21e185SBarry Smith   mg->cyclesperpcapply = n;
13568cc2d5dfSBarry Smith   PetscFunctionReturn(0);
13578cc2d5dfSBarry Smith }
13588cc2d5dfSBarry Smith 
13592134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1360fb15c04eSBarry Smith {
1361fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1362fb15c04eSBarry Smith 
1363fb15c04eSBarry Smith   PetscFunctionBegin;
13642134b1e4SBarry Smith   mg->galerkin = use;
1365fb15c04eSBarry Smith   PetscFunctionReturn(0);
1366fb15c04eSBarry Smith }
1367fb15c04eSBarry Smith 
1368c2be2410SBarry Smith /*@
136997177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
137082b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1371c2be2410SBarry Smith 
1372ad4df100SBarry Smith    Logically Collective on PC
1373c2be2410SBarry Smith 
1374c2be2410SBarry Smith    Input Parameters:
1375c91913d3SJed Brown +  pc - the multigrid context
13762134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1377c2be2410SBarry Smith 
1378c2be2410SBarry Smith    Options Database Key:
13792134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1380c2be2410SBarry Smith 
1381c2be2410SBarry Smith    Level: intermediate
1382c2be2410SBarry Smith 
138395452b02SPatrick Sanan    Notes:
138495452b02SPatrick Sanan     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
13851c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
13861c1aac46SBarry Smith 
13872134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
13883fc8bf9cSBarry Smith 
1389c2be2410SBarry Smith @*/
13902134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1391c2be2410SBarry Smith {
1392fb15c04eSBarry Smith   PetscErrorCode ierr;
1393c2be2410SBarry Smith 
1394c2be2410SBarry Smith   PetscFunctionBegin;
13950700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
13962134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1397c2be2410SBarry Smith   PetscFunctionReturn(0);
1398c2be2410SBarry Smith }
1399c2be2410SBarry Smith 
14003fc8bf9cSBarry Smith /*@
14013fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
140282b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
14033fc8bf9cSBarry Smith 
14043fc8bf9cSBarry Smith    Not Collective
14053fc8bf9cSBarry Smith 
14063fc8bf9cSBarry Smith    Input Parameter:
14073fc8bf9cSBarry Smith .  pc - the multigrid context
14083fc8bf9cSBarry Smith 
14093fc8bf9cSBarry Smith    Output Parameter:
14102134b1e4SBarry 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
14113fc8bf9cSBarry Smith 
14123fc8bf9cSBarry Smith    Level: intermediate
14133fc8bf9cSBarry Smith 
14142134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
14153fc8bf9cSBarry Smith 
14163fc8bf9cSBarry Smith @*/
14172134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
14183fc8bf9cSBarry Smith {
1419f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
14203fc8bf9cSBarry Smith 
14213fc8bf9cSBarry Smith   PetscFunctionBegin;
14220700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
14232134b1e4SBarry Smith   *galerkin = mg->galerkin;
14243fc8bf9cSBarry Smith   PetscFunctionReturn(0);
14253fc8bf9cSBarry Smith }
14263fc8bf9cSBarry Smith 
1427f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1428f3b08a26SMatthew G. Knepley {
1429f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
1430f3b08a26SMatthew G. Knepley 
1431f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1432f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = adapt;
1433f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1434f3b08a26SMatthew G. Knepley }
1435f3b08a26SMatthew G. Knepley 
1436f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1437f3b08a26SMatthew G. Knepley {
1438f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
1439f3b08a26SMatthew G. Knepley 
1440f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1441f3b08a26SMatthew G. Knepley   *adapt = mg->adaptInterpolation;
1442f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1443f3b08a26SMatthew G. Knepley }
1444f3b08a26SMatthew G. Knepley 
144541b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr)
144641b6fd38SMatthew G. Knepley {
144741b6fd38SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
144841b6fd38SMatthew G. Knepley 
144941b6fd38SMatthew G. Knepley   PetscFunctionBegin;
145041b6fd38SMatthew G. Knepley   mg->compatibleRelaxation = cr;
145141b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
145241b6fd38SMatthew G. Knepley }
145341b6fd38SMatthew G. Knepley 
145441b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr)
145541b6fd38SMatthew G. Knepley {
145641b6fd38SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
145741b6fd38SMatthew G. Knepley 
145841b6fd38SMatthew G. Knepley   PetscFunctionBegin;
145941b6fd38SMatthew G. Knepley   *cr = mg->compatibleRelaxation;
146041b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
146141b6fd38SMatthew G. Knepley }
146241b6fd38SMatthew G. Knepley 
1463f3b08a26SMatthew G. Knepley /*@
1464f3b08a26SMatthew 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.
1465f3b08a26SMatthew G. Knepley 
1466f3b08a26SMatthew G. Knepley   Logically Collective on PC
1467f3b08a26SMatthew G. Knepley 
1468f3b08a26SMatthew G. Knepley   Input Parameters:
1469f3b08a26SMatthew G. Knepley + pc    - the multigrid context
1470f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator
1471f3b08a26SMatthew G. Knepley 
1472f3b08a26SMatthew G. Knepley   Options Database Keys:
1473f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp                     - Turn on adaptation
1474f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int>             - The number of modes to use, should be divisible by dimension
1475f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector
1476f3b08a26SMatthew G. Knepley 
1477f3b08a26SMatthew G. Knepley   Level: intermediate
1478f3b08a26SMatthew G. Knepley 
1479f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin
1480f3b08a26SMatthew G. Knepley .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1481f3b08a26SMatthew G. Knepley @*/
1482f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1483f3b08a26SMatthew G. Knepley {
1484f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1485f3b08a26SMatthew G. Knepley 
1486f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1487f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1488f3b08a26SMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));CHKERRQ(ierr);
1489f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1490f3b08a26SMatthew G. Knepley }
1491f3b08a26SMatthew G. Knepley 
1492f3b08a26SMatthew G. Knepley /*@
1493f3b08a26SMatthew 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.
1494f3b08a26SMatthew G. Knepley 
1495f3b08a26SMatthew G. Knepley   Not collective
1496f3b08a26SMatthew G. Knepley 
1497f3b08a26SMatthew G. Knepley   Input Parameter:
1498f3b08a26SMatthew G. Knepley . pc    - the multigrid context
1499f3b08a26SMatthew G. Knepley 
1500f3b08a26SMatthew G. Knepley   Output Parameter:
1501f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator
1502f3b08a26SMatthew G. Knepley 
1503f3b08a26SMatthew G. Knepley   Level: intermediate
1504f3b08a26SMatthew G. Knepley 
1505f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin
1506f3b08a26SMatthew G. Knepley .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1507f3b08a26SMatthew G. Knepley @*/
1508f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1509f3b08a26SMatthew G. Knepley {
1510f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1511f3b08a26SMatthew G. Knepley 
1512f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1513f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1514f3b08a26SMatthew G. Knepley   PetscValidPointer(adapt, 2);
1515f3b08a26SMatthew G. Knepley   ierr = PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));CHKERRQ(ierr);
1516f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1517f3b08a26SMatthew G. Knepley }
1518f3b08a26SMatthew G. Knepley 
15194b9ad928SBarry Smith /*@
152041b6fd38SMatthew G. Knepley   PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation.
152141b6fd38SMatthew G. Knepley 
152241b6fd38SMatthew G. Knepley   Logically Collective on PC
152341b6fd38SMatthew G. Knepley 
152441b6fd38SMatthew G. Knepley   Input Parameters:
152541b6fd38SMatthew G. Knepley + pc - the multigrid context
152641b6fd38SMatthew G. Knepley - cr - flag for compatible relaxation
152741b6fd38SMatthew G. Knepley 
152841b6fd38SMatthew G. Knepley   Options Database Keys:
152941b6fd38SMatthew G. Knepley . -pc_mg_adapt_cr - Turn on compatible relaxation
153041b6fd38SMatthew G. Knepley 
153141b6fd38SMatthew G. Knepley   Level: intermediate
153241b6fd38SMatthew G. Knepley 
153341b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin
153441b6fd38SMatthew G. Knepley .seealso: PCMGGetAdaptCR(), PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
153541b6fd38SMatthew G. Knepley @*/
153641b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr)
153741b6fd38SMatthew G. Knepley {
153841b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
153941b6fd38SMatthew G. Knepley 
154041b6fd38SMatthew G. Knepley   PetscFunctionBegin;
154141b6fd38SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
154241b6fd38SMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr));CHKERRQ(ierr);
154341b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
154441b6fd38SMatthew G. Knepley }
154541b6fd38SMatthew G. Knepley 
154641b6fd38SMatthew G. Knepley /*@
154741b6fd38SMatthew G. Knepley   PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation.
154841b6fd38SMatthew G. Knepley 
154941b6fd38SMatthew G. Knepley   Not collective
155041b6fd38SMatthew G. Knepley 
155141b6fd38SMatthew G. Knepley   Input Parameter:
155241b6fd38SMatthew G. Knepley . pc    - the multigrid context
155341b6fd38SMatthew G. Knepley 
155441b6fd38SMatthew G. Knepley   Output Parameter:
155541b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion
155641b6fd38SMatthew G. Knepley 
155741b6fd38SMatthew G. Knepley   Level: intermediate
155841b6fd38SMatthew G. Knepley 
155941b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin
156041b6fd38SMatthew G. Knepley .seealso: PCMGSetAdaptCR(), PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
156141b6fd38SMatthew G. Knepley @*/
156241b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr)
156341b6fd38SMatthew G. Knepley {
156441b6fd38SMatthew G. Knepley   PetscErrorCode ierr;
156541b6fd38SMatthew G. Knepley 
156641b6fd38SMatthew G. Knepley   PetscFunctionBegin;
156741b6fd38SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
156841b6fd38SMatthew G. Knepley   PetscValidPointer(cr, 2);
156941b6fd38SMatthew G. Knepley   ierr = PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr));CHKERRQ(ierr);
157041b6fd38SMatthew G. Knepley   PetscFunctionReturn(0);
157141b6fd38SMatthew G. Knepley }
157241b6fd38SMatthew G. Knepley 
157341b6fd38SMatthew G. Knepley /*@
157406792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1575710315b6SLawrence Mitchell    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1576710315b6SLawrence Mitchell    pre- and post-smoothing steps.
157706792cafSBarry Smith 
157806792cafSBarry Smith    Logically Collective on PC
157906792cafSBarry Smith 
158006792cafSBarry Smith    Input Parameters:
158106792cafSBarry Smith +  mg - the multigrid context
158206792cafSBarry Smith -  n - the number of smoothing steps
158306792cafSBarry Smith 
158406792cafSBarry Smith    Options Database Key:
1585a2b725a8SWilliam Gropp .  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
158606792cafSBarry Smith 
158706792cafSBarry Smith    Level: advanced
158806792cafSBarry Smith 
158995452b02SPatrick Sanan    Notes:
159095452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
159106792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
159206792cafSBarry Smith 
1593710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp()
159406792cafSBarry Smith @*/
159506792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
159606792cafSBarry Smith {
159706792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
159806792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
159906792cafSBarry Smith   PetscErrorCode ierr;
160006792cafSBarry Smith   PetscInt       i,levels;
160106792cafSBarry Smith 
160206792cafSBarry Smith   PetscFunctionBegin;
160306792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
160406792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
16051a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
160606792cafSBarry Smith   levels = mglevels[0]->levels;
160706792cafSBarry Smith 
160806792cafSBarry Smith   for (i=1; i<levels; i++) {
160906792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
161006792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
161106792cafSBarry Smith     mg->default_smoothu = n;
161206792cafSBarry Smith     mg->default_smoothd = n;
161306792cafSBarry Smith   }
161406792cafSBarry Smith   PetscFunctionReturn(0);
161506792cafSBarry Smith }
161606792cafSBarry Smith 
1617f442ab6aSBarry Smith /*@
16188e5aa403SBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels
1619710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1620f442ab6aSBarry Smith 
1621f442ab6aSBarry Smith    Logically Collective on PC
1622f442ab6aSBarry Smith 
1623f442ab6aSBarry Smith    Input Parameters:
1624f442ab6aSBarry Smith .  pc - the preconditioner context
1625f442ab6aSBarry Smith 
1626f442ab6aSBarry Smith    Options Database Key:
1627f442ab6aSBarry Smith .  -pc_mg_distinct_smoothup
1628f442ab6aSBarry Smith 
1629f442ab6aSBarry Smith    Level: advanced
1630f442ab6aSBarry Smith 
163195452b02SPatrick Sanan    Notes:
163295452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
1633f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1634f442ab6aSBarry Smith 
1635710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth()
1636f442ab6aSBarry Smith @*/
1637f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1638f442ab6aSBarry Smith {
1639f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1640f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1641f442ab6aSBarry Smith   PetscErrorCode ierr;
1642f442ab6aSBarry Smith   PetscInt       i,levels;
1643f442ab6aSBarry Smith   KSP            subksp;
1644f442ab6aSBarry Smith 
1645f442ab6aSBarry Smith   PetscFunctionBegin;
1646f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1647f442ab6aSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1648f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1649f442ab6aSBarry Smith 
1650f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1651710315b6SLawrence Mitchell     const char *prefix = NULL;
1652f442ab6aSBarry Smith     /* make sure smoother up and down are different */
1653f442ab6aSBarry Smith     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1654710315b6SLawrence Mitchell     ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr);
1655710315b6SLawrence Mitchell     ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr);
1656f442ab6aSBarry Smith     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1657f442ab6aSBarry Smith   }
1658f442ab6aSBarry Smith   PetscFunctionReturn(0);
1659f442ab6aSBarry Smith }
1660f442ab6aSBarry Smith 
166107a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
166207a4832bSFande Kong PetscErrorCode  PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[])
166307a4832bSFande Kong {
166407a4832bSFande Kong   PC_MG          *mg        = (PC_MG*)pc->data;
166507a4832bSFande Kong   PC_MG_Levels   **mglevels = mg->levels;
166607a4832bSFande Kong   Mat            *mat;
166707a4832bSFande Kong   PetscInt       l;
166807a4832bSFande Kong   PetscErrorCode ierr;
166907a4832bSFande Kong 
167007a4832bSFande Kong   PetscFunctionBegin;
167107a4832bSFande Kong   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
167207a4832bSFande Kong   ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr);
167307a4832bSFande Kong   for (l=1; l< mg->nlevels; l++) {
167407a4832bSFande Kong     mat[l-1] = mglevels[l]->interpolate;
167507a4832bSFande Kong     ierr = PetscObjectReference((PetscObject)mat[l-1]);CHKERRQ(ierr);
167607a4832bSFande Kong   }
167707a4832bSFande Kong   *num_levels = mg->nlevels;
167807a4832bSFande Kong   *interpolations = mat;
167907a4832bSFande Kong   PetscFunctionReturn(0);
168007a4832bSFande Kong }
168107a4832bSFande Kong 
168207a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
168307a4832bSFande Kong PetscErrorCode  PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[])
168407a4832bSFande Kong {
168507a4832bSFande Kong   PC_MG          *mg        = (PC_MG*)pc->data;
168607a4832bSFande Kong   PC_MG_Levels   **mglevels = mg->levels;
168707a4832bSFande Kong   PetscInt       l;
168807a4832bSFande Kong   Mat            *mat;
168907a4832bSFande Kong   PetscErrorCode ierr;
169007a4832bSFande Kong 
169107a4832bSFande Kong   PetscFunctionBegin;
169207a4832bSFande Kong   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
169307a4832bSFande Kong   ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr);
169407a4832bSFande Kong   for (l=0; l<mg->nlevels-1; l++) {
169507a4832bSFande Kong     ierr = KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));CHKERRQ(ierr);
169607a4832bSFande Kong     ierr = PetscObjectReference((PetscObject)mat[l]);CHKERRQ(ierr);
169707a4832bSFande Kong   }
169807a4832bSFande Kong   *num_levels = mg->nlevels;
169907a4832bSFande Kong   *coarseOperators = mat;
170007a4832bSFande Kong   PetscFunctionReturn(0);
170107a4832bSFande Kong }
170207a4832bSFande Kong 
1703f3b08a26SMatthew G. Knepley /*@C
1704f3b08a26SMatthew G. Knepley   PCMGRegisterCoarseSpaceConstructor -  Adds a method to the PCMG package for coarse space construction.
1705f3b08a26SMatthew G. Knepley 
1706f3b08a26SMatthew G. Knepley   Not collective
1707f3b08a26SMatthew G. Knepley 
1708f3b08a26SMatthew G. Knepley   Input Parameters:
1709f3b08a26SMatthew G. Knepley + name     - name of the constructor
1710f3b08a26SMatthew G. Knepley - function - constructor routine
1711f3b08a26SMatthew G. Knepley 
1712f3b08a26SMatthew G. Knepley   Notes:
1713f3b08a26SMatthew G. Knepley   Calling sequence for the routine:
1714f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1715f3b08a26SMatthew G. Knepley $   pc        - The PC object
1716f3b08a26SMatthew G. Knepley $   l         - The multigrid level, 0 is the coarse level
1717f3b08a26SMatthew G. Knepley $   dm        - The DM for this level
1718f3b08a26SMatthew G. Knepley $   smooth    - The level smoother
1719f3b08a26SMatthew G. Knepley $   Nc        - The size of the coarse space
1720f3b08a26SMatthew G. Knepley $   initGuess - Basis for an initial guess for the space
1721f3b08a26SMatthew G. Knepley $   coarseSp  - A basis for the computed coarse space
1722f3b08a26SMatthew G. Knepley 
1723f3b08a26SMatthew G. Knepley   Level: advanced
1724f3b08a26SMatthew G. Knepley 
1725f3b08a26SMatthew G. Knepley .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister()
1726f3b08a26SMatthew G. Knepley @*/
1727f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1728f3b08a26SMatthew G. Knepley {
1729f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1730f3b08a26SMatthew G. Knepley 
1731f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1732f3b08a26SMatthew G. Knepley   ierr = PCInitializePackage();CHKERRQ(ierr);
1733f3b08a26SMatthew G. Knepley   ierr = PetscFunctionListAdd(&PCMGCoarseList,name,function);CHKERRQ(ierr);
1734f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1735f3b08a26SMatthew G. Knepley }
1736f3b08a26SMatthew G. Knepley 
1737f3b08a26SMatthew G. Knepley /*@C
1738f3b08a26SMatthew G. Knepley   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.
1739f3b08a26SMatthew G. Knepley 
1740f3b08a26SMatthew G. Knepley   Not collective
1741f3b08a26SMatthew G. Knepley 
1742f3b08a26SMatthew G. Knepley   Input Parameter:
1743f3b08a26SMatthew G. Knepley . name     - name of the constructor
1744f3b08a26SMatthew G. Knepley 
1745f3b08a26SMatthew G. Knepley   Output Parameter:
1746f3b08a26SMatthew G. Knepley . function - constructor routine
1747f3b08a26SMatthew G. Knepley 
1748f3b08a26SMatthew G. Knepley   Notes:
1749f3b08a26SMatthew G. Knepley   Calling sequence for the routine:
1750f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1751f3b08a26SMatthew G. Knepley $   pc        - The PC object
1752f3b08a26SMatthew G. Knepley $   l         - The multigrid level, 0 is the coarse level
1753f3b08a26SMatthew G. Knepley $   dm        - The DM for this level
1754f3b08a26SMatthew G. Knepley $   smooth    - The level smoother
1755f3b08a26SMatthew G. Knepley $   Nc        - The size of the coarse space
1756f3b08a26SMatthew G. Knepley $   initGuess - Basis for an initial guess for the space
1757f3b08a26SMatthew G. Knepley $   coarseSp  - A basis for the computed coarse space
1758f3b08a26SMatthew G. Knepley 
1759f3b08a26SMatthew G. Knepley   Level: advanced
1760f3b08a26SMatthew G. Knepley 
1761f3b08a26SMatthew G. Knepley .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister()
1762f3b08a26SMatthew G. Knepley @*/
1763f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1764f3b08a26SMatthew G. Knepley {
1765f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1766f3b08a26SMatthew G. Knepley 
1767f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1768f3b08a26SMatthew G. Knepley   ierr = PetscFunctionListFind(PCMGCoarseList,name,function);CHKERRQ(ierr);
1769f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1770f3b08a26SMatthew G. Knepley }
1771f3b08a26SMatthew G. Knepley 
17724b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
17734b9ad928SBarry Smith 
17743b09bd56SBarry Smith /*MC
1775ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
17763b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
17773b09bd56SBarry Smith 
17783b09bd56SBarry Smith    Options Database Keys:
17793b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1780391689abSStefano Zampini .  -pc_mg_cycle_type <v,w> - provide the cycle desired
17818c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
17823b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1783710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
17842134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
17858cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
17868cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1787e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
17888cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
17898cf6037eSBarry Smith                         to the binary output file called binaryoutput
17903b09bd56SBarry Smith 
179195452b02SPatrick Sanan    Notes:
17920b3c753eSRichard 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
17933b09bd56SBarry Smith 
17948cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
17958cf6037eSBarry Smith 
179623067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
179723067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
179823067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
179923067569SBarry Smith        residual is computed at the end of each cycle.
180023067569SBarry Smith 
18013b09bd56SBarry Smith    Level: intermediate
18023b09bd56SBarry Smith 
18038cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1804710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1805710315b6SLawrence Mitchell            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
180697177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
18070d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
18083b09bd56SBarry Smith M*/
18093b09bd56SBarry Smith 
18108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
18114b9ad928SBarry Smith {
1812f3fbd535SBarry Smith   PC_MG          *mg;
1813f3fbd535SBarry Smith   PetscErrorCode ierr;
1814f3fbd535SBarry Smith 
18154b9ad928SBarry Smith   PetscFunctionBegin;
1816b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1817f3fbd535SBarry Smith   pc->data     = (void*)mg;
1818f3fbd535SBarry Smith   mg->nlevels  = -1;
181910eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
18202134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1821f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = PETSC_FALSE;
1822f3b08a26SMatthew G. Knepley   mg->Nc                 = -1;
1823f3b08a26SMatthew G. Knepley   mg->eigenvalue         = -1;
1824f3fbd535SBarry Smith 
182537a44384SMark Adams   pc->useAmat = PETSC_TRUE;
182637a44384SMark Adams 
18274b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
1828*fcb023d4SJed Brown   pc->ops->applytranspose = PCApplyTranspose_MG;
18294b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1830a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
18314b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
18324b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
18334b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1834fb15c04eSBarry Smith 
1835f3b08a26SMatthew G. Knepley   ierr = PetscObjectComposedDataRegister(&mg->eigenvalue);CHKERRQ(ierr);
1836fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
183797d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
183897d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
1839fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);CHKERRQ(ierr);
1840fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);CHKERRQ(ierr);
1841f3b08a26SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);CHKERRQ(ierr);
1842f3b08a26SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);CHKERRQ(ierr);
184341b6fd38SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG);CHKERRQ(ierr);
184441b6fd38SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG);CHKERRQ(ierr);
18454b9ad928SBarry Smith   PetscFunctionReturn(0);
18464b9ad928SBarry Smith }
1847