xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 391689abd96e70b352f6871f0325c5d36d86735e)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>                    /*I "petscksp.h" I*/
61e25c274SJed Brown #include <petscdm.h>
7*391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*);
84b9ad928SBarry Smith 
931567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
104b9ad928SBarry Smith {
1131567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1231567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
136849ba73SBarry Smith   PetscErrorCode ierr;
1431567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
15a0808db4SHong Zhang   PC             subpc;
16a0808db4SHong Zhang   PCFailedReason pcreason;
174b9ad928SBarry Smith 
184b9ad928SBarry Smith   PetscFunctionBegin;
1963e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2031567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
21a0808db4SHong Zhang   ierr = KSPGetPC(mglevels->smoothd,&subpc);CHKERRQ(ierr);
22a0808db4SHong Zhang   ierr = PCGetSetUpFailedReason(subpc,&pcreason);CHKERRQ(ierr);
23a0808db4SHong Zhang   if (pcreason) {
24a0808db4SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
25a0808db4SHong Zhang   }
2663e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2731567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2863e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2931567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
3063e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
314b9ad928SBarry Smith 
324b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
3331567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
344b9ad928SBarry Smith       PetscReal rnorm;
3531567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
364b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3770441072SBarry Smith         if (rnorm < mg->abstol) {
384d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
3957622a8eSBarry Smith           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr);
404b9ad928SBarry Smith         } else {
414d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
4257622a8eSBarry Smith           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);CHKERRQ(ierr);
434b9ad928SBarry Smith         }
444b9ad928SBarry Smith         PetscFunctionReturn(0);
454b9ad928SBarry Smith       }
464b9ad928SBarry Smith     }
474b9ad928SBarry Smith 
4831567311SBarry Smith     mgc = *(mglevelsin - 1);
4963e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5031567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
5163e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
52efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
534b9ad928SBarry Smith     while (cycles--) {
5431567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
554b9ad928SBarry Smith     }
5663e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5731567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5863e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5963e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
6031567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
6163e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
624b9ad928SBarry Smith   }
634b9ad928SBarry Smith   PetscFunctionReturn(0);
644b9ad928SBarry Smith }
654b9ad928SBarry Smith 
66ace3abfcSBarry 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)
674b9ad928SBarry Smith {
68f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
69f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
70dfbe8321SBarry Smith   PetscErrorCode ierr;
71*391689abSStefano Zampini   PC             tpc;
72*391689abSStefano Zampini   PetscBool      changeu,changed;
73f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
744b9ad928SBarry Smith 
754b9ad928SBarry Smith   PetscFunctionBegin;
76a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
77a762d673SBarry Smith   for (i=0; i<levels; i++) {
78a762d673SBarry Smith     if (!mglevels[i]->A) {
79a762d673SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
80a762d673SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
81a762d673SBarry Smith     }
82a762d673SBarry Smith   }
83*391689abSStefano Zampini 
84*391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
85*391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
86*391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
87*391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
88*391689abSStefano Zampini   if (!changed && !changeu) {
89*391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
90f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
91*391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
92*391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
93*391689abSStefano Zampini       Vec *vec;
94*391689abSStefano Zampini 
95*391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
96*391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
97*391689abSStefano Zampini     }
98*391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
99*391689abSStefano Zampini   }
100f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
1014b9ad928SBarry Smith 
10231567311SBarry Smith   mg->rtol   = rtol;
10331567311SBarry Smith   mg->abstol = abstol;
10431567311SBarry Smith   mg->dtol   = dtol;
1054b9ad928SBarry Smith   if (rtol) {
1064b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
1074b9ad928SBarry Smith     PetscReal rnorm;
1087319c654SBarry Smith     if (zeroguess) {
1097319c654SBarry Smith       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
1107319c654SBarry Smith     } else {
111f3fbd535SBarry Smith       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
1124b9ad928SBarry Smith       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
1137319c654SBarry Smith     }
11431567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
1152fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
1162fa5cd67SKarl Rupp   else mg->ttol = 0.0;
1174b9ad928SBarry Smith 
1184d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
119336babb1SJed Brown      stop prematurely due to small residual */
1204d0a8057SBarry Smith   for (i=1; i<levels; i++) {
121f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
122f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
12323067569SBarry Smith       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
12423067569SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
125f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1264b9ad928SBarry Smith     }
1274d0a8057SBarry Smith   }
1284d0a8057SBarry Smith 
1294d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1304d0a8057SBarry Smith   for (i=0; i<its; i++) {
131f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1324d0a8057SBarry Smith     if (*reason) break;
1334d0a8057SBarry Smith   }
1344d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1354d0a8057SBarry Smith   *outits = i;
136*391689abSStefano Zampini   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
1374b9ad928SBarry Smith   PetscFunctionReturn(0);
1384b9ad928SBarry Smith }
1394b9ad928SBarry Smith 
1403aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1413aeaf226SBarry Smith {
1423aeaf226SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1433aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1443aeaf226SBarry Smith   PetscErrorCode ierr;
1453aeaf226SBarry Smith   PetscInt       i,n;
1463aeaf226SBarry Smith 
1473aeaf226SBarry Smith   PetscFunctionBegin;
1483aeaf226SBarry Smith   if (mglevels) {
1493aeaf226SBarry Smith     n = mglevels[0]->levels;
1503aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1513aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1523aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1533aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
1543aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1553aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
1560de8493bSLawrence Mitchell       ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr);
15773250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1583aeaf226SBarry Smith     }
159*391689abSStefano Zampini     /* this is not null only if the smoother on the finest level
160*391689abSStefano Zampini        changes the rhs during PreSolve */
161*391689abSStefano Zampini     ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr);
1623aeaf226SBarry Smith 
1633aeaf226SBarry Smith     for (i=0; i<n; i++) {
1643aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1653aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1663aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1673aeaf226SBarry Smith       }
1683aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1693aeaf226SBarry Smith     }
1703aeaf226SBarry Smith   }
1713aeaf226SBarry Smith   PetscFunctionReturn(0);
1723aeaf226SBarry Smith }
1733aeaf226SBarry Smith 
1744b9ad928SBarry Smith /*@C
17597177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1764b9ad928SBarry Smith    Must be called before any other MG routine.
1774b9ad928SBarry Smith 
178ad4df100SBarry Smith    Logically Collective on PC
1794b9ad928SBarry Smith 
1804b9ad928SBarry Smith    Input Parameters:
1814b9ad928SBarry Smith +  pc - the preconditioner context
1824b9ad928SBarry Smith .  levels - the number of levels
1834b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1842bf68e3eSBarry Smith            on smaller sets of processors.
1854b9ad928SBarry Smith 
1864b9ad928SBarry Smith    Level: intermediate
1874b9ad928SBarry Smith 
1884b9ad928SBarry Smith    Notes:
1894b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1904b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1914b9ad928SBarry Smith 
1924b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1934b9ad928SBarry Smith 
19497177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1954b9ad928SBarry Smith @*/
1967087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1974b9ad928SBarry Smith {
198dfbe8321SBarry Smith   PetscErrorCode ierr;
199f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
200ce94432eSBarry Smith   MPI_Comm       comm;
2013aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
20210eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
20310167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
204f3fbd535SBarry Smith   PetscInt       i;
205f3fbd535SBarry Smith   PetscMPIInt    size;
206f3fbd535SBarry Smith   const char     *prefix;
207f3fbd535SBarry Smith   PC             ipc;
2083aeaf226SBarry Smith   PetscInt       n;
2094b9ad928SBarry Smith 
2104b9ad928SBarry Smith   PetscFunctionBegin;
2110700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
212548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
213548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
2141a97d7b6SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
2153aeaf226SBarry Smith   if (mglevels) {
21610eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
2173aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
2183aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
2193aeaf226SBarry Smith     n    = mglevels[0]->levels;
2203aeaf226SBarry Smith     for (i=0; i<n; i++) {
2213aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2223aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
2233aeaf226SBarry Smith       }
2243aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
2253aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
2263aeaf226SBarry Smith     }
2273aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
2283aeaf226SBarry Smith   }
229f3fbd535SBarry Smith 
230f3fbd535SBarry Smith   mg->nlevels = levels;
231f3fbd535SBarry Smith 
232785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
2333bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
234f3fbd535SBarry Smith 
235f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
236f3fbd535SBarry Smith 
237a9db87a2SMatthew G Knepley   mg->stageApply = 0;
238f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
239b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
2402fa5cd67SKarl Rupp 
241f3fbd535SBarry Smith     mglevels[i]->level               = i;
242f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
24310eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
244336babb1SJed Brown     mg->default_smoothu              = 2;
245336babb1SJed Brown     mg->default_smoothd              = 2;
24663e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
24763e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
24863e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
24963e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
250f3fbd535SBarry Smith 
251f3fbd535SBarry Smith     if (comms) comm = comms[i];
252f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
253422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
2540ee9a668SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
2550ee9a668SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
2560ee9a668SBarry Smith     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
2570ee9a668SBarry Smith     if (i || levels == 1) {
2580ee9a668SBarry Smith       char tprefix[128];
2590ee9a668SBarry Smith 
260336babb1SJed Brown       ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2610059c7bdSJed Brown       ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
262336babb1SJed Brown       ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
263336babb1SJed Brown       ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
264336babb1SJed Brown       ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
2650ee9a668SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
266f3fbd535SBarry Smith 
2670ee9a668SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
2680ee9a668SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
2690ee9a668SBarry Smith     } else {
270f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
271f3fbd535SBarry Smith 
2729dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
273f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
274f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
275f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
276f3fbd535SBarry Smith       if (size > 1) {
277f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
278f3fbd535SBarry Smith       } else {
279f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
280f3fbd535SBarry Smith       }
281753b7fb9SBarry Smith       ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
282f3fbd535SBarry Smith     }
2833bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
2842fa5cd67SKarl Rupp 
285f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
28631567311SBarry Smith     mg->rtol             = 0.0;
28731567311SBarry Smith     mg->abstol           = 0.0;
28831567311SBarry Smith     mg->dtol             = 0.0;
28931567311SBarry Smith     mg->ttol             = 0.0;
29031567311SBarry Smith     mg->cyclesperpcapply = 1;
291f3fbd535SBarry Smith   }
292f3fbd535SBarry Smith   mg->levels = mglevels;
29310eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2944b9ad928SBarry Smith   PetscFunctionReturn(0);
2954b9ad928SBarry Smith }
2964b9ad928SBarry Smith 
297c07bf074SBarry Smith 
298c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
299c07bf074SBarry Smith {
300c07bf074SBarry Smith   PetscErrorCode ierr;
301a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
302a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
303a06653b4SBarry Smith   PetscInt       i,n;
304c07bf074SBarry Smith 
305c07bf074SBarry Smith   PetscFunctionBegin;
306a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
307a06653b4SBarry Smith   if (mglevels) {
308a06653b4SBarry Smith     n = mglevels[0]->levels;
309a06653b4SBarry Smith     for (i=0; i<n; i++) {
310a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
3116bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
312a06653b4SBarry Smith       }
3136bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
314a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
315a06653b4SBarry Smith     }
316a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
317a06653b4SBarry Smith   }
318c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
319f3fbd535SBarry Smith   PetscFunctionReturn(0);
320f3fbd535SBarry Smith }
321f3fbd535SBarry Smith 
322f3fbd535SBarry Smith 
323f3fbd535SBarry Smith 
32409573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
32509573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
32609573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
327f3fbd535SBarry Smith 
328f3fbd535SBarry Smith /*
329f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
330f3fbd535SBarry Smith              or full cycle of multigrid.
331f3fbd535SBarry Smith 
332f3fbd535SBarry Smith   Note:
333f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
334f3fbd535SBarry Smith */
335f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
336f3fbd535SBarry Smith {
337f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
338f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
339f3fbd535SBarry Smith   PetscErrorCode ierr;
340*391689abSStefano Zampini   PC             tpc;
341f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
342*391689abSStefano Zampini   PetscBool      changeu,changed;
343f3fbd535SBarry Smith 
344f3fbd535SBarry Smith   PetscFunctionBegin;
345a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
346e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
347298cc208SBarry Smith   for (i=0; i<levels; i++) {
348e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
34923ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
350298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
351e1d8e5deSBarry Smith     }
352e1d8e5deSBarry Smith   }
353e1d8e5deSBarry Smith 
354*391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
355*391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
356*391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
357*391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
358*391689abSStefano Zampini   if (!changeu && !changed) {
359*391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
360f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
361*391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
362*391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
363*391689abSStefano Zampini       Vec *vec;
364*391689abSStefano Zampini 
365*391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
366*391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
367*391689abSStefano Zampini     }
368*391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
369*391689abSStefano Zampini   }
370f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
371*391689abSStefano Zampini 
37231567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
373f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
37431567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
3750298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
376f3fbd535SBarry Smith     }
3772fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
37831567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
3792fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
38031567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
3812fa5cd67SKarl Rupp   } else {
382f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
383f3fbd535SBarry Smith   }
384a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
385*391689abSStefano Zampini   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
386f3fbd535SBarry Smith   PetscFunctionReturn(0);
387f3fbd535SBarry Smith }
388f3fbd535SBarry Smith 
389f3fbd535SBarry Smith 
3904416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
391f3fbd535SBarry Smith {
392f3fbd535SBarry Smith   PetscErrorCode   ierr;
393710315b6SLawrence Mitchell   PetscInt         levels,cycles;
3942134b1e4SBarry Smith   PetscBool        flg;
395f3fbd535SBarry Smith   PC_MG            *mg = (PC_MG*)pc->data;
3963d3eaba7SBarry Smith   PC_MG_Levels     **mglevels;
397f3fbd535SBarry Smith   PCMGType         mgtype;
398f3fbd535SBarry Smith   PCMGCycleType    mgctype;
3992134b1e4SBarry Smith   PCMGGalerkinType gtype;
400f3fbd535SBarry Smith 
401f3fbd535SBarry Smith   PetscFunctionBegin;
4021d743356SStefano Zampini   levels = PetscMax(mg->nlevels,1);
403e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
404f3fbd535SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
4051a97d7b6SStefano Zampini   if (!flg && !mg->levels && pc->dm) {
406eb3f98d2SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
407eb3f98d2SBarry Smith     levels++;
4083aeaf226SBarry Smith     mg->usedmfornumberoflevels = PETSC_TRUE;
409eb3f98d2SBarry Smith   }
4100298fd71SBarry Smith   ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
4113aeaf226SBarry Smith   mglevels = mg->levels;
4123aeaf226SBarry Smith 
413f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
414f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
415f3fbd535SBarry Smith   if (flg) {
416f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
4172fa5cd67SKarl Rupp   }
4182134b1e4SBarry Smith   gtype = mg->galerkin;
4192134b1e4SBarry Smith   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
4202134b1e4SBarry Smith   if (flg) {
4212134b1e4SBarry Smith     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
422f3fbd535SBarry Smith   }
423f56b1265SBarry Smith   flg = PETSC_FALSE;
424f442ab6aSBarry Smith   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
425f442ab6aSBarry Smith   if (flg) {
426f442ab6aSBarry Smith     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
427f442ab6aSBarry Smith   }
42831567311SBarry Smith   mgtype = mg->am;
429f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
430f3fbd535SBarry Smith   if (flg) {
431f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
432f3fbd535SBarry Smith   }
43331567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
4345363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
435f3fbd535SBarry Smith     if (flg) {
436f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
437f3fbd535SBarry Smith     }
438f3fbd535SBarry Smith   }
439f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4400298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
441f3fbd535SBarry Smith   if (flg) {
442f3fbd535SBarry Smith     PetscInt i;
443f3fbd535SBarry Smith     char     eventname[128];
4441a97d7b6SStefano Zampini 
445f3fbd535SBarry Smith     levels = mglevels[0]->levels;
446f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
447f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
44863e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
449f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
45063e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
451f3fbd535SBarry Smith       if (i) {
452f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
45363e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
454f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
45563e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
456f3fbd535SBarry Smith       }
457f3fbd535SBarry Smith     }
458eec35531SMatthew G Knepley 
459ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
460eec35531SMatthew G Knepley     {
461eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
462eec35531SMatthew G Knepley       PetscStageLog stageLog;
463eec35531SMatthew G Knepley       PetscInt      st;
464eec35531SMatthew G Knepley 
465eec35531SMatthew G Knepley       PetscFunctionBegin;
466eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
467eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
468eec35531SMatthew G Knepley         PetscBool same;
469eec35531SMatthew G Knepley 
470eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4712fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
472eec35531SMatthew G Knepley       }
473eec35531SMatthew G Knepley       if (!mg->stageApply) {
474eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
475eec35531SMatthew G Knepley       }
476eec35531SMatthew G Knepley     }
477ec60ca73SSatish Balay #endif
478f3fbd535SBarry Smith   }
479f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
480f3fbd535SBarry Smith   PetscFunctionReturn(0);
481f3fbd535SBarry Smith }
482f3fbd535SBarry Smith 
4836a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4846a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
4852134b1e4SBarry Smith const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0};
486f3fbd535SBarry Smith 
4879804daf3SBarry Smith #include <petscdraw.h>
488f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
489f3fbd535SBarry Smith {
490f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
491f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
492f3fbd535SBarry Smith   PetscErrorCode ierr;
493e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
494219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
495f3fbd535SBarry Smith 
496f3fbd535SBarry Smith   PetscFunctionBegin;
497251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4985b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
499219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
500f3fbd535SBarry Smith   if (iascii) {
501e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
502efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
50331567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
50431567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
505f3fbd535SBarry Smith     }
5062134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
507f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
5082134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
5092134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
5102134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
5112134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
5122134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
5132134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
5144f66f45eSBarry Smith     } else {
5154f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
516f3fbd535SBarry Smith     }
5175adeb434SBarry Smith     if (mg->view){
5185adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
5195adeb434SBarry Smith     }
520f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
521f3fbd535SBarry Smith       if (!i) {
52289cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
523f3fbd535SBarry Smith       } else {
52489cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
525f3fbd535SBarry Smith       }
526f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
527f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
528f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
529f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
530f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
531f3fbd535SBarry Smith       } else if (i) {
53289cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
533f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
534f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
535f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
536f3fbd535SBarry Smith       }
537f3fbd535SBarry Smith     }
5385b0b0462SBarry Smith   } else if (isbinary) {
5395b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5405b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5415b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5425b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5435b0b0462SBarry Smith       }
5445b0b0462SBarry Smith     }
545219b1045SBarry Smith   } else if (isdraw) {
546219b1045SBarry Smith     PetscDraw draw;
547b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
548219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
549219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5500298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
551219b1045SBarry Smith     bottom = y - th;
552219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
553b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
554219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
555219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
556219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
557b4375e8dSPeter Brune       } else {
558b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
559b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
560b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
561b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
562b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
563b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
564b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
565b4375e8dSPeter Brune       }
5660298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5671cd381d6SBarry Smith       bottom -= th;
568219b1045SBarry Smith     }
5695b0b0462SBarry Smith   }
570f3fbd535SBarry Smith   PetscFunctionReturn(0);
571f3fbd535SBarry Smith }
572f3fbd535SBarry Smith 
573af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
574af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
575cab2e9ccSBarry Smith 
576f3fbd535SBarry Smith /*
577f3fbd535SBarry Smith     Calls setup for the KSP on each level
578f3fbd535SBarry Smith */
579f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
580f3fbd535SBarry Smith {
581f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
582f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
583f3fbd535SBarry Smith   PetscErrorCode ierr;
5841a97d7b6SStefano Zampini   PetscInt       i,n;
58598e04984SBarry Smith   PC             cpc;
5863db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
587f3fbd535SBarry Smith   Mat            dA,dB;
588f3fbd535SBarry Smith   Vec            tvec;
589218a07d4SBarry Smith   DM             *dms;
590649052a6SBarry Smith   PetscViewer    viewer = 0;
5912134b1e4SBarry Smith   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
592f3fbd535SBarry Smith 
593f3fbd535SBarry Smith   PetscFunctionBegin;
5941a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
5951a97d7b6SStefano Zampini   n = mglevels[0]->levels;
59601bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5973aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5983aeaf226SBarry Smith     PetscInt levels;
5993aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
6003aeaf226SBarry Smith     levels++;
6013aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
6020298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
6033aeaf226SBarry Smith       n        = levels;
6043aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
6053aeaf226SBarry Smith       mglevels =  mg->levels;
6063aeaf226SBarry Smith     }
6073aeaf226SBarry Smith   }
60898e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
6093aeaf226SBarry Smith 
610f3fbd535SBarry Smith 
611f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
612f3fbd535SBarry Smith   /* so use those from global PC */
613f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
6140298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
61598e04984SBarry Smith   if (opsset) {
61698e04984SBarry Smith     Mat mmat;
61723ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
61898e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
61998e04984SBarry Smith   }
6204a5f07a5SMark F. Adams 
6214a5f07a5SMark F. Adams   if (!opsset) {
62271b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
6234a5f07a5SMark F. Adams     if(use_amat){
624f3fbd535SBarry 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);
62523ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
626f3fbd535SBarry Smith     }
6274a5f07a5SMark F. Adams     else {
6284a5f07a5SMark 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);
62923ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
6304a5f07a5SMark F. Adams     }
6314a5f07a5SMark F. Adams   }
632f3fbd535SBarry Smith 
6339c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
6349c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
6359c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
6369c7ffb39SBarry Smith       continue;
6379c7ffb39SBarry Smith     }
6389c7ffb39SBarry Smith   }
6392134b1e4SBarry Smith 
6402134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
6412134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
6422134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
6432134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
6442134b1e4SBarry Smith   }
6452134b1e4SBarry Smith 
6462134b1e4SBarry Smith 
6479c7ffb39SBarry Smith   /*
6489c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6499c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6509c7ffb39SBarry Smith   */
6512134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
6522d2b81a6SBarry Smith 	/* construct the interpolation from the DMs */
6532e499ae9SBarry Smith     Mat p;
65473250ac0SBarry Smith     Vec rscale;
655785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
656218a07d4SBarry Smith     dms[n-1] = pc->dm;
657ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
658ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
65941fce8fdSHong Zhang 	/*
66041fce8fdSHong Zhang 	   Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
66141fce8fdSHong Zhang 	   Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
66241fce8fdSHong Zhang 	   But it is safe to use -dm_mat_type.
66341fce8fdSHong Zhang 
66441fce8fdSHong Zhang 	   The mat type should not be hardcoded like this, we need to find a better way.
66541fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
66641fce8fdSHong Zhang     */
667218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
668942e3340SBarry Smith       DMKSP     kdm;
669eab52d2bSLawrence Mitchell       PetscBool dmhasrestrict, dmhasinject;
6703c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6712134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
672942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
673d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
674d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6750298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6760298fd71SBarry Smith       kdm->rhsctx          = NULL;
67724c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
678e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6792d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
68000ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
68173250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6826bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
683218a07d4SBarry Smith       }
6843ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
6853ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
6863ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
6873ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
6883ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
6893ad4599aSBarry Smith       }
690eab52d2bSLawrence Mitchell       ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr);
691eab52d2bSLawrence Mitchell       if (dmhasinject && !mglevels[i+1]->inject){
692eab52d2bSLawrence Mitchell         ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr);
693eab52d2bSLawrence Mitchell         ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr);
694eab52d2bSLawrence Mitchell         ierr = MatDestroy(&p);CHKERRQ(ierr);
695eab52d2bSLawrence Mitchell       }
69624c3aa18SBarry Smith     }
6972d2b81a6SBarry Smith 
698ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
6992d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
700ce4cda84SJed Brown   }
701cab2e9ccSBarry Smith 
702ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
7032134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
704cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
705cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
706218a07d4SBarry Smith   }
707218a07d4SBarry Smith 
7082134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
7092134b1e4SBarry Smith     Mat       A,B;
7102134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
7112134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
7122134b1e4SBarry Smith 
7132134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
7142134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
7152134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
716f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
717b9367914SBarry 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");
718b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
719b9367914SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
720b9367914SBarry Smith       }
721b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
722b9367914SBarry Smith         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
723b9367914SBarry Smith       }
7242134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
7252134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
7262134b1e4SBarry Smith       }
7272134b1e4SBarry Smith       if (doA) {
7282df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
7292134b1e4SBarry Smith       }
7302134b1e4SBarry Smith       if (doB) {
7312df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
732b9367914SBarry Smith       }
7332134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
7342134b1e4SBarry Smith       if (!doA && dAeqdB) {
7352134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
7362134b1e4SBarry Smith         A = B;
7372134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX ) {
7382134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
7392134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
740b9367914SBarry Smith       }
7412134b1e4SBarry Smith       if (!doB && dAeqdB) {
7422134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
7432134b1e4SBarry Smith         B = A;
7442134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
74523ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
7462134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
7477d537d92SJed Brown       }
7482134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
7492134b1e4SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
7502134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
7512134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
7522134b1e4SBarry Smith       }
7532134b1e4SBarry Smith       dA = A;
754f3fbd535SBarry Smith       dB = B;
755f3fbd535SBarry Smith     }
756f3fbd535SBarry Smith   }
7572134b1e4SBarry Smith   if (needRestricts && pc->dm && pc->dm->x) {
758cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
759cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
760c88c5224SJed Brown       Mat R;
761c88c5224SJed Brown       Vec rscale;
762cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
763cab2e9ccSBarry Smith         Vec *vecs;
7642a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
765cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
766cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
767cab2e9ccSBarry Smith       }
768c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
769c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
770c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
771c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
772cab2e9ccSBarry Smith     }
773f3fbd535SBarry Smith   }
7742134b1e4SBarry Smith   if (needRestricts && pc->dm) {
775caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
776ccceb50cSJed Brown       DM  dmfine,dmcoarse;
777caa4e7f2SJed Brown       Mat Restrict,Inject;
778caa4e7f2SJed Brown       Vec rscale;
779ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
780ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
781caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
782caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
783eab52d2bSLawrence Mitchell       ierr   = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr);
784caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
785caa4e7f2SJed Brown     }
786caa4e7f2SJed Brown   }
787f3fbd535SBarry Smith 
788f3fbd535SBarry Smith   if (!pc->setupcalled) {
789f3fbd535SBarry Smith     for (i=0; i<n; i++) {
790f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
791f3fbd535SBarry Smith     }
792f3fbd535SBarry Smith     for (i=1; i<n; i++) {
793f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
794f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
795f3fbd535SBarry Smith       }
796f3fbd535SBarry Smith     }
7973ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
798f3fbd535SBarry Smith     for (i=1; i<n; i++) {
7993ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
8003ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
801f3fbd535SBarry Smith     }
802f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
803f3fbd535SBarry Smith       if (!mglevels[i]->b) {
804f3fbd535SBarry Smith         Vec *vec;
8052a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
806f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
8076bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
808f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
809f3fbd535SBarry Smith       }
810f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
811f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
812f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
8136bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
814f3fbd535SBarry Smith       }
815f3fbd535SBarry Smith       if (!mglevels[i]->x) {
816f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
817f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
8186bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
819f3fbd535SBarry Smith       }
820f3fbd535SBarry Smith     }
821f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
822f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
823f3fbd535SBarry Smith       Vec *vec;
8242a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
825f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
8266bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
827f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
828f3fbd535SBarry Smith     }
829f3fbd535SBarry Smith   }
830f3fbd535SBarry Smith 
83198e04984SBarry Smith   if (pc->dm) {
83298e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
83398e04984SBarry Smith     for (i=0; i<n-1; i++) {
83498e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
83598e04984SBarry Smith     }
83698e04984SBarry Smith   }
837f3fbd535SBarry Smith 
838f3fbd535SBarry Smith   for (i=1; i<n; i++) {
8392a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
840f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
841f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
842f3fbd535SBarry Smith     }
84363e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
844f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
845899639b0SHong Zhang     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
846899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
847899639b0SHong Zhang     }
84863e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
849d42688cbSBarry Smith     if (!mglevels[i]->residual) {
850d42688cbSBarry Smith       Mat mat;
85113842ffbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
85254b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
853d42688cbSBarry Smith     }
854f3fbd535SBarry Smith   }
855f3fbd535SBarry Smith   for (i=1; i<n; i++) {
856f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
857f3fbd535SBarry Smith       Mat downmat,downpmat;
858f3fbd535SBarry Smith 
859f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
8600298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
861f3fbd535SBarry Smith       if (!opsset) {
86223ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
86323ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
864f3fbd535SBarry Smith       }
865f3fbd535SBarry Smith 
866f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
86763e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
868f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
869899639b0SHong Zhang       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) {
870899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
871899639b0SHong Zhang       }
87263e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
873f3fbd535SBarry Smith     }
874f3fbd535SBarry Smith   }
875f3fbd535SBarry Smith 
87663e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
877f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
878899639b0SHong Zhang   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
879899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
880899639b0SHong Zhang   }
88163e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
882f3fbd535SBarry Smith 
883f3fbd535SBarry Smith   /*
884f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
885e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
886f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
887f3fbd535SBarry Smith 
888f3fbd535SBarry Smith    Only support one or the other at the same time.
889f3fbd535SBarry Smith   */
890f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
891c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
892ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
893f3fbd535SBarry Smith   dump = PETSC_FALSE;
894f3fbd535SBarry Smith #endif
895c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
896ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
897f3fbd535SBarry Smith 
898f3fbd535SBarry Smith   if (viewer) {
899f3fbd535SBarry Smith     for (i=1; i<n; i++) {
900f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
901f3fbd535SBarry Smith     }
902f3fbd535SBarry Smith     for (i=0; i<n; i++) {
903f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
904f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
905f3fbd535SBarry Smith     }
906f3fbd535SBarry Smith   }
907f3fbd535SBarry Smith   PetscFunctionReturn(0);
908f3fbd535SBarry Smith }
909f3fbd535SBarry Smith 
910f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
911f3fbd535SBarry Smith 
9124b9ad928SBarry Smith /*@
91397177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
9144b9ad928SBarry Smith 
9154b9ad928SBarry Smith    Not Collective
9164b9ad928SBarry Smith 
9174b9ad928SBarry Smith    Input Parameter:
9184b9ad928SBarry Smith .  pc - the preconditioner context
9194b9ad928SBarry Smith 
9204b9ad928SBarry Smith    Output parameter:
9214b9ad928SBarry Smith .  levels - the number of levels
9224b9ad928SBarry Smith 
9234b9ad928SBarry Smith    Level: advanced
9244b9ad928SBarry Smith 
9254b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
9264b9ad928SBarry Smith 
92797177400SBarry Smith .seealso: PCMGSetLevels()
9284b9ad928SBarry Smith @*/
9297087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
9304b9ad928SBarry Smith {
931f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9324b9ad928SBarry Smith 
9334b9ad928SBarry Smith   PetscFunctionBegin;
9340700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9354482741eSBarry Smith   PetscValidIntPointer(levels,2);
936f3fbd535SBarry Smith   *levels = mg->nlevels;
9374b9ad928SBarry Smith   PetscFunctionReturn(0);
9384b9ad928SBarry Smith }
9394b9ad928SBarry Smith 
9404b9ad928SBarry Smith /*@
94197177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
9424b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
9434b9ad928SBarry Smith 
944ad4df100SBarry Smith    Logically Collective on PC
9454b9ad928SBarry Smith 
9464b9ad928SBarry Smith    Input Parameters:
9474b9ad928SBarry Smith +  pc - the preconditioner context
9489dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
9499dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
9504b9ad928SBarry Smith 
9514b9ad928SBarry Smith    Options Database Key:
9524b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
9534b9ad928SBarry Smith    additive, full, kaskade
9544b9ad928SBarry Smith 
9554b9ad928SBarry Smith    Level: advanced
9564b9ad928SBarry Smith 
9574b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
9584b9ad928SBarry Smith 
95997177400SBarry Smith .seealso: PCMGSetLevels()
9604b9ad928SBarry Smith @*/
9617087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
9624b9ad928SBarry Smith {
963f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9644b9ad928SBarry Smith 
9654b9ad928SBarry Smith   PetscFunctionBegin;
9660700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
967c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
96831567311SBarry Smith   mg->am = form;
9699dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
970c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
971c60c7ad4SBarry Smith   PetscFunctionReturn(0);
972c60c7ad4SBarry Smith }
973c60c7ad4SBarry Smith 
974c60c7ad4SBarry Smith /*@
975c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
976c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
977c60c7ad4SBarry Smith 
978c60c7ad4SBarry Smith    Logically Collective on PC
979c60c7ad4SBarry Smith 
980c60c7ad4SBarry Smith    Input Parameter:
981c60c7ad4SBarry Smith .  pc - the preconditioner context
982c60c7ad4SBarry Smith 
983c60c7ad4SBarry Smith    Output Parameter:
984c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
985c60c7ad4SBarry Smith 
986c60c7ad4SBarry Smith 
987c60c7ad4SBarry Smith    Level: advanced
988c60c7ad4SBarry Smith 
989c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
990c60c7ad4SBarry Smith 
991c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
992c60c7ad4SBarry Smith @*/
993c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
994c60c7ad4SBarry Smith {
995c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
996c60c7ad4SBarry Smith 
997c60c7ad4SBarry Smith   PetscFunctionBegin;
998c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
999c60c7ad4SBarry Smith   *type = mg->am;
10004b9ad928SBarry Smith   PetscFunctionReturn(0);
10014b9ad928SBarry Smith }
10024b9ad928SBarry Smith 
10034b9ad928SBarry Smith /*@
10040d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
10054b9ad928SBarry Smith    complicated cycling.
10064b9ad928SBarry Smith 
1007ad4df100SBarry Smith    Logically Collective on PC
10084b9ad928SBarry Smith 
10094b9ad928SBarry Smith    Input Parameters:
1010c2be2410SBarry Smith +  pc - the multigrid context
1011c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
10124b9ad928SBarry Smith 
10134b9ad928SBarry Smith    Options Database Key:
1014c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
10154b9ad928SBarry Smith 
10164b9ad928SBarry Smith    Level: advanced
10174b9ad928SBarry Smith 
10184b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10194b9ad928SBarry Smith 
10200d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
10214b9ad928SBarry Smith @*/
10227087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
10234b9ad928SBarry Smith {
1024f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
1025f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
102679416396SBarry Smith   PetscInt     i,levels;
10274b9ad928SBarry Smith 
10284b9ad928SBarry Smith   PetscFunctionBegin;
10290700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
103069fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
10311a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1032f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10332fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
10344b9ad928SBarry Smith   PetscFunctionReturn(0);
10354b9ad928SBarry Smith }
10364b9ad928SBarry Smith 
10378cc2d5dfSBarry Smith /*@
10388cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
10398cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
10408cc2d5dfSBarry Smith 
1041ad4df100SBarry Smith    Logically Collective on PC
10428cc2d5dfSBarry Smith 
10438cc2d5dfSBarry Smith    Input Parameters:
10448cc2d5dfSBarry Smith +  pc - the multigrid context
10458cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
10468cc2d5dfSBarry Smith 
10478cc2d5dfSBarry Smith    Options Database Key:
1048e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
10498cc2d5dfSBarry Smith 
10508cc2d5dfSBarry Smith    Level: advanced
10518cc2d5dfSBarry Smith 
10528cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
10538cc2d5dfSBarry Smith 
10548cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10558cc2d5dfSBarry Smith 
10568cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
10578cc2d5dfSBarry Smith @*/
10587087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
10598cc2d5dfSBarry Smith {
1060f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
10618cc2d5dfSBarry Smith 
10628cc2d5dfSBarry Smith   PetscFunctionBegin;
10630700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1064c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
10652a21e185SBarry Smith   mg->cyclesperpcapply = n;
10668cc2d5dfSBarry Smith   PetscFunctionReturn(0);
10678cc2d5dfSBarry Smith }
10688cc2d5dfSBarry Smith 
10692134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1070fb15c04eSBarry Smith {
1071fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1072fb15c04eSBarry Smith 
1073fb15c04eSBarry Smith   PetscFunctionBegin;
10742134b1e4SBarry Smith   mg->galerkin = use;
1075fb15c04eSBarry Smith   PetscFunctionReturn(0);
1076fb15c04eSBarry Smith }
1077fb15c04eSBarry Smith 
1078c2be2410SBarry Smith /*@
107997177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
108082b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1081c2be2410SBarry Smith 
1082ad4df100SBarry Smith    Logically Collective on PC
1083c2be2410SBarry Smith 
1084c2be2410SBarry Smith    Input Parameters:
1085c91913d3SJed Brown +  pc - the multigrid context
10862134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1087c2be2410SBarry Smith 
1088c2be2410SBarry Smith    Options Database Key:
10892134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1090c2be2410SBarry Smith 
1091c2be2410SBarry Smith    Level: intermediate
1092c2be2410SBarry Smith 
10931c1aac46SBarry Smith    Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
10941c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
10951c1aac46SBarry Smith 
1096c2be2410SBarry Smith .keywords: MG, set, Galerkin
1097c2be2410SBarry Smith 
10982134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
10993fc8bf9cSBarry Smith 
1100c2be2410SBarry Smith @*/
11012134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1102c2be2410SBarry Smith {
1103fb15c04eSBarry Smith   PetscErrorCode ierr;
1104c2be2410SBarry Smith 
1105c2be2410SBarry Smith   PetscFunctionBegin;
11060700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11072134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1108c2be2410SBarry Smith   PetscFunctionReturn(0);
1109c2be2410SBarry Smith }
1110c2be2410SBarry Smith 
11113fc8bf9cSBarry Smith /*@
11123fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
111382b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
11143fc8bf9cSBarry Smith 
11153fc8bf9cSBarry Smith    Not Collective
11163fc8bf9cSBarry Smith 
11173fc8bf9cSBarry Smith    Input Parameter:
11183fc8bf9cSBarry Smith .  pc - the multigrid context
11193fc8bf9cSBarry Smith 
11203fc8bf9cSBarry Smith    Output Parameter:
11212134b1e4SBarry 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
11223fc8bf9cSBarry Smith 
11233fc8bf9cSBarry Smith    Level: intermediate
11243fc8bf9cSBarry Smith 
11253fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
11263fc8bf9cSBarry Smith 
11272134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
11283fc8bf9cSBarry Smith 
11293fc8bf9cSBarry Smith @*/
11302134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
11313fc8bf9cSBarry Smith {
1132f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
11333fc8bf9cSBarry Smith 
11343fc8bf9cSBarry Smith   PetscFunctionBegin;
11350700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11362134b1e4SBarry Smith   *galerkin = mg->galerkin;
11373fc8bf9cSBarry Smith   PetscFunctionReturn(0);
11383fc8bf9cSBarry Smith }
11393fc8bf9cSBarry Smith 
11404b9ad928SBarry Smith /*@
114106792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1142710315b6SLawrence Mitchell    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1143710315b6SLawrence Mitchell    pre- and post-smoothing steps.
114406792cafSBarry Smith 
114506792cafSBarry Smith    Logically Collective on PC
114606792cafSBarry Smith 
114706792cafSBarry Smith    Input Parameters:
114806792cafSBarry Smith +  mg - the multigrid context
114906792cafSBarry Smith -  n - the number of smoothing steps
115006792cafSBarry Smith 
115106792cafSBarry Smith    Options Database Key:
115206792cafSBarry Smith +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
115306792cafSBarry Smith 
115406792cafSBarry Smith    Level: advanced
115506792cafSBarry Smith 
115606792cafSBarry Smith    Notes: this does not set a value on the coarsest grid, since we assume that
115706792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
115806792cafSBarry Smith 
115906792cafSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
116006792cafSBarry Smith 
1161710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp()
116206792cafSBarry Smith @*/
116306792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
116406792cafSBarry Smith {
116506792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
116606792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
116706792cafSBarry Smith   PetscErrorCode ierr;
116806792cafSBarry Smith   PetscInt       i,levels;
116906792cafSBarry Smith 
117006792cafSBarry Smith   PetscFunctionBegin;
117106792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
117206792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
11731a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
117406792cafSBarry Smith   levels = mglevels[0]->levels;
117506792cafSBarry Smith 
117606792cafSBarry Smith   for (i=1; i<levels; i++) {
117706792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
117806792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
117906792cafSBarry Smith     mg->default_smoothu = n;
118006792cafSBarry Smith     mg->default_smoothd = n;
118106792cafSBarry Smith   }
118206792cafSBarry Smith   PetscFunctionReturn(0);
118306792cafSBarry Smith }
118406792cafSBarry Smith 
1185f442ab6aSBarry Smith /*@
1186f442ab6aSBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels
1187710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1188f442ab6aSBarry Smith 
1189f442ab6aSBarry Smith    Logically Collective on PC
1190f442ab6aSBarry Smith 
1191f442ab6aSBarry Smith    Input Parameters:
1192f442ab6aSBarry Smith .  pc - the preconditioner context
1193f442ab6aSBarry Smith 
1194f442ab6aSBarry Smith    Options Database Key:
1195f442ab6aSBarry Smith .  -pc_mg_distinct_smoothup
1196f442ab6aSBarry Smith 
1197f442ab6aSBarry Smith    Level: advanced
1198f442ab6aSBarry Smith 
1199f442ab6aSBarry Smith    Notes: this does not set a value on the coarsest grid, since we assume that
1200f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1201f442ab6aSBarry Smith 
1202f442ab6aSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1203f442ab6aSBarry Smith 
1204710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth()
1205f442ab6aSBarry Smith @*/
1206f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1207f442ab6aSBarry Smith {
1208f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1209f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1210f442ab6aSBarry Smith   PetscErrorCode ierr;
1211f442ab6aSBarry Smith   PetscInt       i,levels;
1212f442ab6aSBarry Smith   KSP            subksp;
1213f442ab6aSBarry Smith 
1214f442ab6aSBarry Smith   PetscFunctionBegin;
1215f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1216f442ab6aSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1217f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1218f442ab6aSBarry Smith 
1219f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1220710315b6SLawrence Mitchell     const char *prefix = NULL;
1221f442ab6aSBarry Smith     /* make sure smoother up and down are different */
1222f442ab6aSBarry Smith     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1223710315b6SLawrence Mitchell     ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr);
1224710315b6SLawrence Mitchell     ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr);
1225f442ab6aSBarry Smith     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1226f442ab6aSBarry Smith   }
1227f442ab6aSBarry Smith   PetscFunctionReturn(0);
1228f442ab6aSBarry Smith }
1229f442ab6aSBarry Smith 
12304b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
12314b9ad928SBarry Smith 
12323b09bd56SBarry Smith /*MC
1233ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
12343b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
12353b09bd56SBarry Smith 
12363b09bd56SBarry Smith    Options Database Keys:
12373b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1238*391689abSStefano Zampini .  -pc_mg_cycle_type <v,w> - provide the cycle desired
12398c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
12403b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1241710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
12422134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
12438cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
12448cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1245e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
12468cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
12478cf6037eSBarry Smith                         to the binary output file called binaryoutput
12483b09bd56SBarry Smith 
1249e941fc33SBarry Smith    Notes: If one uses a Krylov method such GMRES or CG as the smoother than one must use KSPFGMRES, KSPGCG, or KSPRICHARDSON as the outer Krylov method
12503b09bd56SBarry Smith 
12518cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
12528cf6037eSBarry Smith 
125323067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
125423067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
125523067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
125623067569SBarry Smith        residual is computed at the end of each cycle.
125723067569SBarry Smith 
12583b09bd56SBarry Smith    Level: intermediate
12593b09bd56SBarry Smith 
12608f87f92bSBarry Smith    Concepts: multigrid/multilevel
12613b09bd56SBarry Smith 
12628cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1263710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1264710315b6SLawrence Mitchell            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
126597177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12660d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12673b09bd56SBarry Smith M*/
12683b09bd56SBarry Smith 
12698cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
12704b9ad928SBarry Smith {
1271f3fbd535SBarry Smith   PC_MG          *mg;
1272f3fbd535SBarry Smith   PetscErrorCode ierr;
1273f3fbd535SBarry Smith 
12744b9ad928SBarry Smith   PetscFunctionBegin;
1275b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1276f3fbd535SBarry Smith   pc->data     = (void*)mg;
1277f3fbd535SBarry Smith   mg->nlevels  = -1;
127810eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
12792134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1280f3fbd535SBarry Smith 
128137a44384SMark Adams   pc->useAmat = PETSC_TRUE;
128237a44384SMark Adams 
12834b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
12844b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1285a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
12864b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
12874b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
12884b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1289fb15c04eSBarry Smith 
1290fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
12914b9ad928SBarry Smith   PetscFunctionReturn(0);
12924b9ad928SBarry Smith }
1293