xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision c31cb41c0ee9b115779dca6896353b82fd361af9)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
54b9ad928SBarry Smith */
67c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
74b9ad928SBarry Smith 
84b9ad928SBarry Smith 
94b9ad928SBarry Smith #undef __FUNCT__
109dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private"
1131567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
124b9ad928SBarry Smith {
1331567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1431567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
156849ba73SBarry Smith   PetscErrorCode ierr;
1631567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
174b9ad928SBarry Smith 
184b9ad928SBarry Smith   PetscFunctionBegin;
194b9ad928SBarry Smith 
2063e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2131567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
2263e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2331567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2463e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2531567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2663e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
274b9ad928SBarry Smith 
284b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
2931567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
304b9ad928SBarry Smith       PetscReal rnorm;
3131567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
324b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3370441072SBarry Smith         if (rnorm < mg->abstol) {
344d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
351e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
364b9ad928SBarry Smith         } else {
374d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
381e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr);
394b9ad928SBarry Smith         }
404b9ad928SBarry Smith         PetscFunctionReturn(0);
414b9ad928SBarry Smith       }
424b9ad928SBarry Smith     }
434b9ad928SBarry Smith 
4431567311SBarry Smith     mgc = *(mglevelsin - 1);
4563e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4631567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
4763e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
48efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
494b9ad928SBarry Smith     while (cycles--) {
5031567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
514b9ad928SBarry Smith     }
5263e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5331567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5463e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5563e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5631567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
5763e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
584b9ad928SBarry Smith   }
594b9ad928SBarry Smith   PetscFunctionReturn(0);
604b9ad928SBarry Smith }
614b9ad928SBarry Smith 
624b9ad928SBarry Smith #undef __FUNCT__
634b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
64ace3abfcSBarry 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)
654b9ad928SBarry Smith {
66f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
67f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
68dfbe8321SBarry Smith   PetscErrorCode ierr;
69f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
704b9ad928SBarry Smith 
714b9ad928SBarry Smith   PetscFunctionBegin;
72f3fbd535SBarry Smith   mglevels[levels-1]->b    = b;
73f3fbd535SBarry Smith   mglevels[levels-1]->x    = x;
744b9ad928SBarry Smith 
7531567311SBarry Smith   mg->rtol = rtol;
7631567311SBarry Smith   mg->abstol = abstol;
7731567311SBarry Smith   mg->dtol = dtol;
784b9ad928SBarry Smith   if (rtol) {
794b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
804b9ad928SBarry Smith     PetscReal rnorm;
817319c654SBarry Smith     if (zeroguess) {
827319c654SBarry Smith       ierr               = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
837319c654SBarry Smith     } else {
84f3fbd535SBarry Smith       ierr               = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
854b9ad928SBarry Smith       ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
867319c654SBarry Smith     }
8731567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
8870441072SBarry Smith   } else if (abstol) {
8931567311SBarry Smith     mg->ttol = abstol;
904b9ad928SBarry Smith   } else {
9131567311SBarry Smith     mg->ttol = 0.0;
924b9ad928SBarry Smith   }
934b9ad928SBarry Smith 
944d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
954d0a8057SBarry Smith      stop prematurely do to small residual */
964d0a8057SBarry Smith   for (i=1; i<levels; i++) {
97f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
98f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
99f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1004b9ad928SBarry Smith     }
1014d0a8057SBarry Smith   }
1024d0a8057SBarry Smith 
1034d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1044d0a8057SBarry Smith   for (i=0; i<its; i++) {
105f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1064d0a8057SBarry Smith     if (*reason) break;
1074d0a8057SBarry Smith   }
1084d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1094d0a8057SBarry Smith   *outits = i;
1104b9ad928SBarry Smith   PetscFunctionReturn(0);
1114b9ad928SBarry Smith }
1124b9ad928SBarry Smith 
1134b9ad928SBarry Smith #undef __FUNCT__
1149dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
1154b9ad928SBarry Smith /*@C
11697177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1174b9ad928SBarry Smith    Must be called before any other MG routine.
1184b9ad928SBarry Smith 
119ad4df100SBarry Smith    Logically Collective on PC
1204b9ad928SBarry Smith 
1214b9ad928SBarry Smith    Input Parameters:
1224b9ad928SBarry Smith +  pc - the preconditioner context
1234b9ad928SBarry Smith .  levels - the number of levels
1244b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1254b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
1264b9ad928SBarry Smith 
1274b9ad928SBarry Smith    Level: intermediate
1284b9ad928SBarry Smith 
1294b9ad928SBarry Smith    Notes:
1304b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1314b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1324b9ad928SBarry Smith 
1334b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1344b9ad928SBarry Smith 
13597177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1364b9ad928SBarry Smith @*/
1377087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1384b9ad928SBarry Smith {
139dfbe8321SBarry Smith   PetscErrorCode ierr;
140f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
141f3fbd535SBarry Smith   MPI_Comm       comm = ((PetscObject)pc)->comm;
142f3fbd535SBarry Smith   PC_MG_Levels   **mglevels;
143f3fbd535SBarry Smith   PetscInt       i;
144f3fbd535SBarry Smith   PetscMPIInt    size;
145f3fbd535SBarry Smith   const char     *prefix;
146f3fbd535SBarry Smith   PC             ipc;
1474b9ad928SBarry Smith 
1484b9ad928SBarry Smith   PetscFunctionBegin;
1490700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
150e7e72b3dSBarry Smith   if (mg->nlevels > -1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Number levels already set for MG\n  make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
151e32f2f54SBarry Smith   if (mg->levels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc, this array should not yet exist");
152c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,levels,2);
153f3fbd535SBarry Smith 
154f3fbd535SBarry Smith   mg->nlevels      = levels;
155218a07d4SBarry Smith   mg->galerkin     = PETSC_FALSE;
156218a07d4SBarry Smith   mg->galerkinused = PETSC_FALSE;
157f3fbd535SBarry Smith 
158f3fbd535SBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr);
159f3fbd535SBarry Smith   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
160f3fbd535SBarry Smith 
161f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
162f3fbd535SBarry Smith 
163f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
164f3fbd535SBarry Smith     ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr);
165f3fbd535SBarry Smith     mglevels[i]->level           = i;
166f3fbd535SBarry Smith     mglevels[i]->levels          = levels;
167f3fbd535SBarry Smith     mglevels[i]->cycles          = PC_MG_CYCLE_V;
16831567311SBarry Smith     mg->default_smoothu = 1;
16931567311SBarry Smith     mg->default_smoothd = 1;
17063e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
17163e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
17263e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
17363e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
174f3fbd535SBarry Smith 
175f3fbd535SBarry Smith     if (comms) comm = comms[i];
176f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
177f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
17831567311SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
179f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
180f3fbd535SBarry Smith 
181f3fbd535SBarry Smith     /* do special stuff for coarse grid */
182f3fbd535SBarry Smith     if (!i && levels > 1) {
183f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
184f3fbd535SBarry Smith 
185f3fbd535SBarry Smith       /* coarse solve is (redundant) LU by default */
186f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
187f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
188f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
189f3fbd535SBarry Smith       if (size > 1) {
190f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
191f3fbd535SBarry Smith       } else {
192f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
193f3fbd535SBarry Smith       }
194f3fbd535SBarry Smith 
195f3fbd535SBarry Smith     } else {
196f3fbd535SBarry Smith       char tprefix[128];
197f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
198f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
199f3fbd535SBarry Smith     }
200f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr);
201f3fbd535SBarry Smith     mglevels[i]->smoothu    = mglevels[i]->smoothd;
20231567311SBarry Smith     mg->rtol                = 0.0;
20331567311SBarry Smith     mg->abstol              = 0.0;
20431567311SBarry Smith     mg->dtol                = 0.0;
20531567311SBarry Smith     mg->ttol                = 0.0;
20631567311SBarry Smith     mg->cyclesperpcapply    = 1;
207f3fbd535SBarry Smith   }
20831567311SBarry Smith   mg->am          = PC_MG_MULTIPLICATIVE;
209f3fbd535SBarry Smith   mg->levels      = mglevels;
2104b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
2114b9ad928SBarry Smith   PetscFunctionReturn(0);
2124b9ad928SBarry Smith }
2134b9ad928SBarry Smith 
2144b9ad928SBarry Smith #undef __FUNCT__
215c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG_Private"
216c07bf074SBarry Smith PetscErrorCode PCDestroy_MG_Private(PC pc)
217f3fbd535SBarry Smith {
218f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
219f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
220f3fbd535SBarry Smith   PetscErrorCode ierr;
221f3fbd535SBarry Smith   PetscInt       i,n;
222f3fbd535SBarry Smith 
223f3fbd535SBarry Smith   PetscFunctionBegin;
224f3fbd535SBarry Smith   if (mglevels) {
225f3fbd535SBarry Smith     n = mglevels[0]->levels;
226f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
227f3fbd535SBarry Smith       if (mglevels[i+1]->r) {ierr = VecDestroy(mglevels[i+1]->r);CHKERRQ(ierr);}
228f3fbd535SBarry Smith       if (mglevels[i]->b) {ierr = VecDestroy(mglevels[i]->b);CHKERRQ(ierr);}
229f3fbd535SBarry Smith       if (mglevels[i]->x) {ierr = VecDestroy(mglevels[i]->x);CHKERRQ(ierr);}
230f3fbd535SBarry Smith       if (mglevels[i+1]->restrct) {ierr = MatDestroy(mglevels[i+1]->restrct);CHKERRQ(ierr);}
231f3fbd535SBarry Smith       if (mglevels[i+1]->interpolate) {ierr = MatDestroy(mglevels[i+1]->interpolate);CHKERRQ(ierr);}
232f3fbd535SBarry Smith     }
233f3fbd535SBarry Smith 
234f3fbd535SBarry Smith     for (i=0; i<n; i++) {
235f3fbd535SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
236f3fbd535SBarry Smith 	ierr = KSPDestroy(mglevels[i]->smoothd);CHKERRQ(ierr);
237f3fbd535SBarry Smith       }
238f3fbd535SBarry Smith       ierr = KSPDestroy(mglevels[i]->smoothu);CHKERRQ(ierr);
239f3fbd535SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
240f3fbd535SBarry Smith     }
241*c31cb41cSBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
242f3fbd535SBarry Smith   }
243c07bf074SBarry Smith   mg->nlevels = -1;
244c07bf074SBarry Smith   PetscFunctionReturn(0);
245c07bf074SBarry Smith }
246c07bf074SBarry Smith 
247c07bf074SBarry Smith #undef __FUNCT__
248c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
249c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
250c07bf074SBarry Smith {
251c07bf074SBarry Smith   PetscErrorCode ierr;
252c07bf074SBarry Smith 
253c07bf074SBarry Smith   PetscFunctionBegin;
254c07bf074SBarry Smith   ierr = PCDestroy_MG_Private(pc);CHKERRQ(ierr);
255*c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
256f3fbd535SBarry Smith   PetscFunctionReturn(0);
257f3fbd535SBarry Smith }
258f3fbd535SBarry Smith 
259f3fbd535SBarry Smith 
260f3fbd535SBarry Smith 
26109573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
26209573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
26309573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
264f3fbd535SBarry Smith 
265f3fbd535SBarry Smith /*
266f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
267f3fbd535SBarry Smith              or full cycle of multigrid.
268f3fbd535SBarry Smith 
269f3fbd535SBarry Smith   Note:
270f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
271f3fbd535SBarry Smith */
272f3fbd535SBarry Smith #undef __FUNCT__
273f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
274f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
275f3fbd535SBarry Smith {
276f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
277f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
278f3fbd535SBarry Smith   PetscErrorCode ierr;
279f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
280f3fbd535SBarry Smith 
281f3fbd535SBarry Smith   PetscFunctionBegin;
282e1d8e5deSBarry Smith 
283e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
284e1d8e5deSBarry Smith   for (i=0; i<levels-1; i++) {
285e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
286e1d8e5deSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
287e1d8e5deSBarry Smith     }
288e1d8e5deSBarry Smith   }
289e1d8e5deSBarry Smith 
290f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
291f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
29231567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
293f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
29431567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
295f3fbd535SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
296f3fbd535SBarry Smith     }
297f3fbd535SBarry Smith   }
29831567311SBarry Smith   else if (mg->am == PC_MG_ADDITIVE) {
29931567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
300f3fbd535SBarry Smith   }
30131567311SBarry Smith   else if (mg->am == PC_MG_KASKADE) {
30231567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
303f3fbd535SBarry Smith   }
304f3fbd535SBarry Smith   else {
305f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
306f3fbd535SBarry Smith   }
307f3fbd535SBarry Smith   PetscFunctionReturn(0);
308f3fbd535SBarry Smith }
309f3fbd535SBarry Smith 
310f3fbd535SBarry Smith 
311f3fbd535SBarry Smith #undef __FUNCT__
312f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
313f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc)
314f3fbd535SBarry Smith {
315f3fbd535SBarry Smith   PetscErrorCode ierr;
316f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
317ace3abfcSBarry Smith   PetscBool      flg;
318f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
319f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
320f3fbd535SBarry Smith   PCMGType       mgtype;
321f3fbd535SBarry Smith   PCMGCycleType  mgctype;
322f3fbd535SBarry Smith 
323f3fbd535SBarry Smith   PetscFunctionBegin;
324f3fbd535SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
32518aabeadSBarry Smith     if (!mglevels) {
326f3fbd535SBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
327f3fbd535SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
328f3fbd535SBarry Smith       mglevels = mg->levels;
329f3fbd535SBarry Smith     }
330f3fbd535SBarry Smith     mgctype = (PCMGCycleType) mglevels[0]->cycles;
331f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
332f3fbd535SBarry Smith     if (flg) {
333f3fbd535SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
334f3fbd535SBarry Smith     };
335f3fbd535SBarry Smith     flg  = PETSC_FALSE;
336acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
337f3fbd535SBarry Smith     if (flg) {
338f3fbd535SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
339f3fbd535SBarry Smith     }
340f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
341f3fbd535SBarry Smith     if (flg) {
342f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
343f3fbd535SBarry Smith     }
344f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
345f3fbd535SBarry Smith     if (flg) {
346f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
347f3fbd535SBarry Smith     }
34831567311SBarry Smith     mgtype = mg->am;
349f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
350f3fbd535SBarry Smith     if (flg) {
351f3fbd535SBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
352f3fbd535SBarry Smith     }
35331567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
35431567311SBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
355f3fbd535SBarry Smith       if (flg) {
356f3fbd535SBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
357f3fbd535SBarry Smith       }
358f3fbd535SBarry Smith     }
359f3fbd535SBarry Smith     flg  = PETSC_FALSE;
360acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
361f3fbd535SBarry Smith     if (flg) {
362f3fbd535SBarry Smith       PetscInt i;
363f3fbd535SBarry Smith       char     eventname[128];
36463e6d426SJed Brown       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
365f3fbd535SBarry Smith       levels = mglevels[0]->levels;
366f3fbd535SBarry Smith       for (i=0; i<levels; i++) {
367f3fbd535SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
36863e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
369f3fbd535SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
37063e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
371f3fbd535SBarry Smith         if (i) {
372f3fbd535SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
37363e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
374f3fbd535SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
37563e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
376f3fbd535SBarry Smith         }
377f3fbd535SBarry Smith       }
378f3fbd535SBarry Smith     }
379f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
380f3fbd535SBarry Smith   PetscFunctionReturn(0);
381f3fbd535SBarry Smith }
382f3fbd535SBarry Smith 
383f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
384f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
385f3fbd535SBarry Smith 
386f3fbd535SBarry Smith #undef __FUNCT__
387f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
388f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
389f3fbd535SBarry Smith {
390f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
391f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
392f3fbd535SBarry Smith   PetscErrorCode ierr;
393f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
394ace3abfcSBarry Smith   PetscBool      iascii;
395f3fbd535SBarry Smith 
396f3fbd535SBarry Smith   PetscFunctionBegin;
3972692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
398f3fbd535SBarry Smith   if (iascii) {
39931567311SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,(mglevels[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr);
40031567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
40131567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
402f3fbd535SBarry Smith     }
403218a07d4SBarry Smith     if (mg->galerkin) {
404f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4054f66f45eSBarry Smith     } else {
4064f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
407f3fbd535SBarry Smith     }
408f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
409f3fbd535SBarry Smith       if (!i) {
4102d19a89fSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
411f3fbd535SBarry Smith       } else {
41231567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
413f3fbd535SBarry Smith       }
414f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
415f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
416f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
417f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
418f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
419f3fbd535SBarry Smith       } else if (i){
42031567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr);
421f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
422f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
423f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
424f3fbd535SBarry Smith       }
425f3fbd535SBarry Smith     }
426f3fbd535SBarry Smith   } else {
42765e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
428f3fbd535SBarry Smith   }
429f3fbd535SBarry Smith   PetscFunctionReturn(0);
430f3fbd535SBarry Smith }
431f3fbd535SBarry Smith 
432f3fbd535SBarry Smith /*
433f3fbd535SBarry Smith     Calls setup for the KSP on each level
434f3fbd535SBarry Smith */
435f3fbd535SBarry Smith #undef __FUNCT__
436f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
437f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
438f3fbd535SBarry Smith {
439f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
440f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
441f3fbd535SBarry Smith   PetscErrorCode          ierr;
442f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
443f3fbd535SBarry Smith   PC                      cpc,mpc;
444ace3abfcSBarry Smith   PetscBool               preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset;
445f3fbd535SBarry Smith   PetscViewerASCIIMonitor ascii;
446f3fbd535SBarry Smith   PetscViewer             viewer = PETSC_NULL;
447f3fbd535SBarry Smith   MPI_Comm                comm;
448f3fbd535SBarry Smith   Mat                     dA,dB;
449f3fbd535SBarry Smith   MatStructure            uflag;
450f3fbd535SBarry Smith   Vec                     tvec;
451218a07d4SBarry Smith   DM                      *dms;
452f3fbd535SBarry Smith 
453f3fbd535SBarry Smith   PetscFunctionBegin;
454f3fbd535SBarry Smith 
455f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
456f3fbd535SBarry Smith   /* so use those from global PC */
457f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
458f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
459f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
460f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr);
4611338a6b9SJed Brown   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) {
462f3fbd535SBarry 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);
463f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
464f3fbd535SBarry Smith   }
465f3fbd535SBarry Smith 
4662d2b81a6SBarry Smith   if (pc->dm && !pc->setupcalled) {
4672d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
4682e499ae9SBarry Smith     Mat p;
469218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
470218a07d4SBarry Smith     dms[n-1] = pc->dm;
471218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
472218a07d4SBarry Smith       ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr);
47311629dbeSBarry Smith       ierr = DMSetFunction(dms[i],0);
47411629dbeSBarry Smith       ierr = DMSetInitialGuess(dms[i],0);
47524c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
4762d2b81a6SBarry Smith 	ierr = DMGetInterpolation(dms[i],dms[i+1],&p,PETSC_NULL);CHKERRQ(ierr);
4772d2b81a6SBarry Smith 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
4782d2b81a6SBarry Smith         ierr = MatDestroy(p);CHKERRQ(ierr);
479218a07d4SBarry Smith       }
48024c3aa18SBarry Smith     }
4812d2b81a6SBarry Smith 
4822d2b81a6SBarry Smith     if (!mg->galerkin) {
4832e499ae9SBarry Smith       /* each coarse level gets its DM */
4842d2b81a6SBarry Smith       for (i=n-2; i>-1; i--) {
4852e499ae9SBarry Smith         ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
4862d2b81a6SBarry Smith       }
4872d2b81a6SBarry Smith     }
4882d2b81a6SBarry Smith 
489218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
490218a07d4SBarry Smith       ierr = DMDestroy(dms[i]);CHKERRQ(ierr);
491218a07d4SBarry Smith     }
4922d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
493218a07d4SBarry Smith   }
494218a07d4SBarry Smith 
495218a07d4SBarry Smith   if (mg->galerkin) {
496f3fbd535SBarry Smith     Mat B;
497218a07d4SBarry Smith     mg->galerkinused = PETSC_TRUE;
498f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
499f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
500f3fbd535SBarry Smith     if (!pc->setupcalled) {
501f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
502f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
503f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
504f3fbd535SBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
505f3fbd535SBarry Smith         dB   = B;
506f3fbd535SBarry Smith       }
507cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
508f3fbd535SBarry Smith     } else {
509f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
510f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
511f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
512f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
513f3fbd535SBarry Smith         dB   = B;
514f3fbd535SBarry Smith       }
515f3fbd535SBarry Smith     }
516f3fbd535SBarry Smith   }
517f3fbd535SBarry Smith 
518f3fbd535SBarry Smith   if (!pc->setupcalled) {
519acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr);
520f3fbd535SBarry Smith 
521f3fbd535SBarry Smith     for (i=0; i<n; i++) {
522f3fbd535SBarry Smith       if (monitor) {
523f3fbd535SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr);
524f3fbd535SBarry Smith         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
525f3fbd535SBarry Smith         ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
526f3fbd535SBarry Smith       }
527f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
528f3fbd535SBarry Smith     }
529f3fbd535SBarry Smith     for (i=1; i<n; i++) {
530f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
531f3fbd535SBarry Smith         if (monitor) {
532f3fbd535SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr);
533f3fbd535SBarry Smith           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
534f3fbd535SBarry Smith           ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
535f3fbd535SBarry Smith         }
536f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
537f3fbd535SBarry Smith       }
538f3fbd535SBarry Smith     }
539f3fbd535SBarry Smith     for (i=1; i<n; i++) {
540f3fbd535SBarry Smith       if (!mglevels[i]->residual) {
541f3fbd535SBarry Smith         Mat mat;
542f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
543f3fbd535SBarry Smith         ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
544f3fbd535SBarry Smith       }
545f3fbd535SBarry Smith       if (mglevels[i]->restrct && !mglevels[i]->interpolate) {
546f3fbd535SBarry Smith         ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr);
547f3fbd535SBarry Smith       }
548f3fbd535SBarry Smith       if (!mglevels[i]->restrct && mglevels[i]->interpolate) {
549f3fbd535SBarry Smith         ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr);
550f3fbd535SBarry Smith       }
551f3fbd535SBarry Smith #if defined(PETSC_USE_DEBUG)
552f3fbd535SBarry Smith       if (!mglevels[i]->restrct || !mglevels[i]->interpolate) {
55365e19b50SBarry Smith         SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
554f3fbd535SBarry Smith       }
555f3fbd535SBarry Smith #endif
556f3fbd535SBarry Smith     }
557f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
558f3fbd535SBarry Smith       if (!mglevels[i]->b) {
559f3fbd535SBarry Smith         Vec *vec;
560f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
561f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
562f3fbd535SBarry Smith         ierr = VecDestroy(*vec);CHKERRQ(ierr);
563f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
564f3fbd535SBarry Smith       }
565f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
566f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
567f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
568f3fbd535SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
569f3fbd535SBarry Smith       }
570f3fbd535SBarry Smith       if (!mglevels[i]->x) {
571f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
572f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
573f3fbd535SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
574f3fbd535SBarry Smith       }
575f3fbd535SBarry Smith     }
576f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
577f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
578f3fbd535SBarry Smith       Vec *vec;
579f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
580f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
581f3fbd535SBarry Smith       ierr = VecDestroy(*vec);CHKERRQ(ierr);
582f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
583f3fbd535SBarry Smith     }
584f3fbd535SBarry Smith   }
585f3fbd535SBarry Smith 
586f3fbd535SBarry Smith 
587f3fbd535SBarry Smith   for (i=1; i<n; i++) {
588f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
589f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
590f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
591f3fbd535SBarry Smith     }
59263e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
593f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
59463e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
595f3fbd535SBarry Smith   }
596f3fbd535SBarry Smith   for (i=1; i<n; i++) {
597f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
598f3fbd535SBarry Smith       Mat          downmat,downpmat;
599f3fbd535SBarry Smith       MatStructure matflag;
600ace3abfcSBarry Smith       PetscBool    opsset;
601f3fbd535SBarry Smith 
602f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
603f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
604f3fbd535SBarry Smith       if (!opsset) {
605f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
606f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
607f3fbd535SBarry Smith       }
608f3fbd535SBarry Smith 
609f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
61063e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
611f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
61263e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
613f3fbd535SBarry Smith     }
614f3fbd535SBarry Smith   }
615f3fbd535SBarry Smith 
616f3fbd535SBarry Smith   /*
617f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
618f3fbd535SBarry Smith   */
619f3fbd535SBarry Smith   ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
620f3fbd535SBarry Smith   if (preonly) {
621f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
622f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
623f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
624f3fbd535SBarry Smith     if (!lu && !redundant && !cholesky) {
625f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
626f3fbd535SBarry Smith     }
627f3fbd535SBarry Smith   }
628f3fbd535SBarry Smith 
629f3fbd535SBarry Smith   if (!pc->setupcalled) {
630f3fbd535SBarry Smith     if (monitor) {
631f3fbd535SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr);
632f3fbd535SBarry Smith       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
633f3fbd535SBarry Smith       ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
634f3fbd535SBarry Smith     }
635f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
636f3fbd535SBarry Smith   }
637f3fbd535SBarry Smith 
63863e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
639f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
64063e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
641f3fbd535SBarry Smith 
642f3fbd535SBarry Smith   /*
643f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
644e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
645f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
646f3fbd535SBarry Smith 
647f3fbd535SBarry Smith    Only support one or the other at the same time.
648f3fbd535SBarry Smith   */
649f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
650acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
651f3fbd535SBarry Smith   if (dump) {
652f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
653f3fbd535SBarry Smith   }
654f3fbd535SBarry Smith   dump = PETSC_FALSE;
655f3fbd535SBarry Smith #endif
656acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
657f3fbd535SBarry Smith   if (dump) {
658f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
659f3fbd535SBarry Smith   }
660f3fbd535SBarry Smith 
661f3fbd535SBarry Smith   if (viewer) {
662f3fbd535SBarry Smith     for (i=1; i<n; i++) {
663f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
664f3fbd535SBarry Smith     }
665f3fbd535SBarry Smith     for (i=0; i<n; i++) {
666f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
667f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
668f3fbd535SBarry Smith     }
669f3fbd535SBarry Smith   }
670f3fbd535SBarry Smith   PetscFunctionReturn(0);
671f3fbd535SBarry Smith }
672f3fbd535SBarry Smith 
673f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
674f3fbd535SBarry Smith 
675f3fbd535SBarry Smith #undef __FUNCT__
6769dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
6774b9ad928SBarry Smith /*@
67897177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
6794b9ad928SBarry Smith 
6804b9ad928SBarry Smith    Not Collective
6814b9ad928SBarry Smith 
6824b9ad928SBarry Smith    Input Parameter:
6834b9ad928SBarry Smith .  pc - the preconditioner context
6844b9ad928SBarry Smith 
6854b9ad928SBarry Smith    Output parameter:
6864b9ad928SBarry Smith .  levels - the number of levels
6874b9ad928SBarry Smith 
6884b9ad928SBarry Smith    Level: advanced
6894b9ad928SBarry Smith 
6904b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
6914b9ad928SBarry Smith 
69297177400SBarry Smith .seealso: PCMGSetLevels()
6934b9ad928SBarry Smith @*/
6947087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
6954b9ad928SBarry Smith {
696f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
6974b9ad928SBarry Smith 
6984b9ad928SBarry Smith   PetscFunctionBegin;
6990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
7004482741eSBarry Smith   PetscValidIntPointer(levels,2);
701f3fbd535SBarry Smith   *levels = mg->nlevels;
7024b9ad928SBarry Smith   PetscFunctionReturn(0);
7034b9ad928SBarry Smith }
7044b9ad928SBarry Smith 
7054b9ad928SBarry Smith #undef __FUNCT__
7069dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
7074b9ad928SBarry Smith /*@
70897177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
7094b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
7104b9ad928SBarry Smith 
711ad4df100SBarry Smith    Logically Collective on PC
7124b9ad928SBarry Smith 
7134b9ad928SBarry Smith    Input Parameters:
7144b9ad928SBarry Smith +  pc - the preconditioner context
7159dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
7169dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
7174b9ad928SBarry Smith 
7184b9ad928SBarry Smith    Options Database Key:
7194b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
7204b9ad928SBarry Smith    additive, full, kaskade
7214b9ad928SBarry Smith 
7224b9ad928SBarry Smith    Level: advanced
7234b9ad928SBarry Smith 
7244b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
7254b9ad928SBarry Smith 
72697177400SBarry Smith .seealso: PCMGSetLevels()
7274b9ad928SBarry Smith @*/
7287087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
7294b9ad928SBarry Smith {
730f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
7314b9ad928SBarry Smith 
7324b9ad928SBarry Smith   PetscFunctionBegin;
7330700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
734c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
73531567311SBarry Smith   mg->am = form;
7369dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
7374b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
7384b9ad928SBarry Smith   PetscFunctionReturn(0);
7394b9ad928SBarry Smith }
7404b9ad928SBarry Smith 
7414b9ad928SBarry Smith #undef __FUNCT__
7420d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
7434b9ad928SBarry Smith /*@
7440d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
7454b9ad928SBarry Smith    complicated cycling.
7464b9ad928SBarry Smith 
747ad4df100SBarry Smith    Logically Collective on PC
7484b9ad928SBarry Smith 
7494b9ad928SBarry Smith    Input Parameters:
750c2be2410SBarry Smith +  pc - the multigrid context
7510d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
7524b9ad928SBarry Smith 
7534b9ad928SBarry Smith    Options Database Key:
7540d353602SBarry Smith $  -pc_mg_cycle_type v or w
7554b9ad928SBarry Smith 
7564b9ad928SBarry Smith    Level: advanced
7574b9ad928SBarry Smith 
7584b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7594b9ad928SBarry Smith 
7600d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
7614b9ad928SBarry Smith @*/
7627087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
7634b9ad928SBarry Smith {
764f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
765f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
76679416396SBarry Smith   PetscInt     i,levels;
7674b9ad928SBarry Smith 
7684b9ad928SBarry Smith   PetscFunctionBegin;
7690700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
770e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
771c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
772f3fbd535SBarry Smith   levels = mglevels[0]->levels;
7734b9ad928SBarry Smith 
7744b9ad928SBarry Smith   for (i=0; i<levels; i++) {
775f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
7764b9ad928SBarry Smith   }
7774b9ad928SBarry Smith   PetscFunctionReturn(0);
7784b9ad928SBarry Smith }
7794b9ad928SBarry Smith 
7804b9ad928SBarry Smith #undef __FUNCT__
7818cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
7828cc2d5dfSBarry Smith /*@
7838cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
7848cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
7858cc2d5dfSBarry Smith 
786ad4df100SBarry Smith    Logically Collective on PC
7878cc2d5dfSBarry Smith 
7888cc2d5dfSBarry Smith    Input Parameters:
7898cc2d5dfSBarry Smith +  pc - the multigrid context
7908cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
7918cc2d5dfSBarry Smith 
7928cc2d5dfSBarry Smith    Options Database Key:
7938cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
7948cc2d5dfSBarry Smith 
7958cc2d5dfSBarry Smith    Level: advanced
7968cc2d5dfSBarry Smith 
7978cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
7988cc2d5dfSBarry Smith 
7998cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8008cc2d5dfSBarry Smith 
8018cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
8028cc2d5dfSBarry Smith @*/
8037087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
8048cc2d5dfSBarry Smith {
805f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
806f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
8078cc2d5dfSBarry Smith   PetscInt     i,levels;
8088cc2d5dfSBarry Smith 
8098cc2d5dfSBarry Smith   PetscFunctionBegin;
8100700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
811e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
812c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
813f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8148cc2d5dfSBarry Smith 
8158cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
81631567311SBarry Smith     mg->cyclesperpcapply  = n;
8178cc2d5dfSBarry Smith   }
8188cc2d5dfSBarry Smith   PetscFunctionReturn(0);
8198cc2d5dfSBarry Smith }
8208cc2d5dfSBarry Smith 
8218cc2d5dfSBarry Smith #undef __FUNCT__
8229dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
823c2be2410SBarry Smith /*@
82497177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
825c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
826c2be2410SBarry Smith 
827ad4df100SBarry Smith    Logically Collective on PC
828c2be2410SBarry Smith 
829c2be2410SBarry Smith    Input Parameters:
8303fc8bf9cSBarry Smith .  pc - the multigrid context
831c2be2410SBarry Smith 
832c2be2410SBarry Smith    Options Database Key:
833c2be2410SBarry Smith $  -pc_mg_galerkin
834c2be2410SBarry Smith 
835c2be2410SBarry Smith    Level: intermediate
836c2be2410SBarry Smith 
837c2be2410SBarry Smith .keywords: MG, set, Galerkin
838c2be2410SBarry Smith 
8393fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
8403fc8bf9cSBarry Smith 
841c2be2410SBarry Smith @*/
8427087cfbeSBarry Smith PetscErrorCode  PCMGSetGalerkin(PC pc)
843c2be2410SBarry Smith {
844f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
845c2be2410SBarry Smith 
846c2be2410SBarry Smith   PetscFunctionBegin;
8470700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
848218a07d4SBarry Smith   mg->galerkin = PETSC_TRUE;
849c2be2410SBarry Smith   PetscFunctionReturn(0);
850c2be2410SBarry Smith }
851c2be2410SBarry Smith 
852c2be2410SBarry Smith #undef __FUNCT__
8533fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
8543fc8bf9cSBarry Smith /*@
8553fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
8563fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
8573fc8bf9cSBarry Smith 
8583fc8bf9cSBarry Smith    Not Collective
8593fc8bf9cSBarry Smith 
8603fc8bf9cSBarry Smith    Input Parameter:
8613fc8bf9cSBarry Smith .  pc - the multigrid context
8623fc8bf9cSBarry Smith 
8633fc8bf9cSBarry Smith    Output Parameter:
8643fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
8653fc8bf9cSBarry Smith 
8663fc8bf9cSBarry Smith    Options Database Key:
8673fc8bf9cSBarry Smith $  -pc_mg_galerkin
8683fc8bf9cSBarry Smith 
8693fc8bf9cSBarry Smith    Level: intermediate
8703fc8bf9cSBarry Smith 
8713fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
8723fc8bf9cSBarry Smith 
8733fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
8743fc8bf9cSBarry Smith 
8753fc8bf9cSBarry Smith @*/
8767087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
8773fc8bf9cSBarry Smith {
878f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
8793fc8bf9cSBarry Smith 
8803fc8bf9cSBarry Smith   PetscFunctionBegin;
8810700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
882218a07d4SBarry Smith   *galerkin = mg->galerkin;
8833fc8bf9cSBarry Smith   PetscFunctionReturn(0);
8843fc8bf9cSBarry Smith }
8853fc8bf9cSBarry Smith 
8863fc8bf9cSBarry Smith #undef __FUNCT__
8879dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
8884b9ad928SBarry Smith /*@
88997177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
89097177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
8914b9ad928SBarry Smith    pre-smoothing steps on different levels.
8924b9ad928SBarry Smith 
893ad4df100SBarry Smith    Logically Collective on PC
8944b9ad928SBarry Smith 
8954b9ad928SBarry Smith    Input Parameters:
8964b9ad928SBarry Smith +  mg - the multigrid context
8974b9ad928SBarry Smith -  n - the number of smoothing steps
8984b9ad928SBarry Smith 
8994b9ad928SBarry Smith    Options Database Key:
9004b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
9014b9ad928SBarry Smith 
9024b9ad928SBarry Smith    Level: advanced
9034b9ad928SBarry Smith 
9044b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
9054b9ad928SBarry Smith 
90697177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
9074b9ad928SBarry Smith @*/
9087087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
9094b9ad928SBarry Smith {
910f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
911f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9126849ba73SBarry Smith   PetscErrorCode ierr;
91379416396SBarry Smith   PetscInt       i,levels;
9144b9ad928SBarry Smith 
9154b9ad928SBarry Smith   PetscFunctionBegin;
9160700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
917e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
918c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
919f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9204b9ad928SBarry Smith 
921b05257ddSBarry Smith   for (i=1; i<levels; i++) {
9224b9ad928SBarry Smith     /* make sure smoother up and down are different */
92397177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
924f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
92531567311SBarry Smith     mg->default_smoothd = n;
9264b9ad928SBarry Smith   }
9274b9ad928SBarry Smith   PetscFunctionReturn(0);
9284b9ad928SBarry Smith }
9294b9ad928SBarry Smith 
9304b9ad928SBarry Smith #undef __FUNCT__
9319dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
9324b9ad928SBarry Smith /*@
93397177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
93497177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
9354b9ad928SBarry Smith    post-smoothing steps on different levels.
9364b9ad928SBarry Smith 
937ad4df100SBarry Smith    Logically Collective on PC
9384b9ad928SBarry Smith 
9394b9ad928SBarry Smith    Input Parameters:
9404b9ad928SBarry Smith +  mg - the multigrid context
9414b9ad928SBarry Smith -  n - the number of smoothing steps
9424b9ad928SBarry Smith 
9434b9ad928SBarry Smith    Options Database Key:
9444b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
9454b9ad928SBarry Smith 
9464b9ad928SBarry Smith    Level: advanced
9474b9ad928SBarry Smith 
9484b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
949a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
9504b9ad928SBarry Smith 
9514b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
9524b9ad928SBarry Smith 
95397177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
9544b9ad928SBarry Smith @*/
9557087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
9564b9ad928SBarry Smith {
957f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
958f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9596849ba73SBarry Smith   PetscErrorCode ierr;
96079416396SBarry Smith   PetscInt       i,levels;
9614b9ad928SBarry Smith 
9624b9ad928SBarry Smith   PetscFunctionBegin;
9630700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
964e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
965c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
966f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9674b9ad928SBarry Smith 
9684b9ad928SBarry Smith   for (i=1; i<levels; i++) {
9694b9ad928SBarry Smith     /* make sure smoother up and down are different */
97097177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
971f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
97231567311SBarry Smith     mg->default_smoothu = n;
9734b9ad928SBarry Smith   }
9744b9ad928SBarry Smith   PetscFunctionReturn(0);
9754b9ad928SBarry Smith }
9764b9ad928SBarry Smith 
9774b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
9784b9ad928SBarry Smith 
9793b09bd56SBarry Smith /*MC
980ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
9813b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
9823b09bd56SBarry Smith 
9833b09bd56SBarry Smith    Options Database Keys:
9843b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
9850d353602SBarry Smith .  -pc_mg_cycles v or w
98679416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
9873b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
9883b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
9893b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
9903b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
99168eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
9923b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
993e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
9943b09bd56SBarry Smith 
99524c3aa18SBarry Smith    Notes: By default this uses GMRES on the fine grid smoother so this should be used with KSPFGMRES or the smoother changed to not use GMRES
9963b09bd56SBarry Smith 
9973b09bd56SBarry Smith    Level: intermediate
9983b09bd56SBarry Smith 
9998f87f92bSBarry Smith    Concepts: multigrid/multilevel
10003b09bd56SBarry Smith 
100124c3aa18SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC
10020d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
100397177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
100497177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
10050d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
10063b09bd56SBarry Smith M*/
10073b09bd56SBarry Smith 
10084b9ad928SBarry Smith EXTERN_C_BEGIN
10094b9ad928SBarry Smith #undef __FUNCT__
10104b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
10117087cfbeSBarry Smith PetscErrorCode  PCCreate_MG(PC pc)
10124b9ad928SBarry Smith {
1013f3fbd535SBarry Smith   PC_MG          *mg;
1014f3fbd535SBarry Smith   PetscErrorCode ierr;
1015f3fbd535SBarry Smith 
10164b9ad928SBarry Smith   PetscFunctionBegin;
1017f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1018f3fbd535SBarry Smith   pc->data    = (void*)mg;
1019f3fbd535SBarry Smith   mg->nlevels = -1;
1020f3fbd535SBarry Smith 
10214b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
10224b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
10234b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
10244b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
10254b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
10264b9ad928SBarry Smith   PetscFunctionReturn(0);
10274b9ad928SBarry Smith }
10284b9ad928SBarry Smith EXTERN_C_END
1029