xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 548767e0e3f20f644ed4e14175e5b52c15b3dbd2)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscpcmg.h" I*/
64b9ad928SBarry Smith 
74b9ad928SBarry Smith 
84b9ad928SBarry Smith #undef __FUNCT__
99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private"
1031567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
114b9ad928SBarry Smith {
1231567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1331567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
146849ba73SBarry Smith   PetscErrorCode ierr;
1531567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
184b9ad928SBarry Smith 
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 */
2163e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2231567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2363e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2431567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2563e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
264b9ad928SBarry Smith 
274b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
2831567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
294b9ad928SBarry Smith       PetscReal rnorm;
3031567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
314b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3270441072SBarry Smith         if (rnorm < mg->abstol) {
334d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
341e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
354b9ad928SBarry Smith         } else {
364d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
371e2582c4SBarry 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);
384b9ad928SBarry Smith         }
394b9ad928SBarry Smith         PetscFunctionReturn(0);
404b9ad928SBarry Smith       }
414b9ad928SBarry Smith     }
424b9ad928SBarry Smith 
4331567311SBarry Smith     mgc = *(mglevelsin - 1);
4463e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4531567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
4663e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
47efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
484b9ad928SBarry Smith     while (cycles--) {
4931567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
504b9ad928SBarry Smith     }
5163e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5231567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5363e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5463e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5531567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
5663e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith   }
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
614b9ad928SBarry Smith #undef __FUNCT__
624b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
63ace3abfcSBarry 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)
644b9ad928SBarry Smith {
65f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
66f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
67dfbe8321SBarry Smith   PetscErrorCode ierr;
68f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
694b9ad928SBarry Smith 
704b9ad928SBarry Smith   PetscFunctionBegin;
71f3fbd535SBarry Smith   mglevels[levels-1]->b    = b;
72f3fbd535SBarry Smith   mglevels[levels-1]->x    = x;
734b9ad928SBarry Smith 
7431567311SBarry Smith   mg->rtol = rtol;
7531567311SBarry Smith   mg->abstol = abstol;
7631567311SBarry Smith   mg->dtol = dtol;
774b9ad928SBarry Smith   if (rtol) {
784b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
794b9ad928SBarry Smith     PetscReal rnorm;
807319c654SBarry Smith     if (zeroguess) {
817319c654SBarry Smith       ierr               = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
827319c654SBarry Smith     } else {
83f3fbd535SBarry Smith       ierr               = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
844b9ad928SBarry Smith       ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
857319c654SBarry Smith     }
8631567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
8770441072SBarry Smith   } else if (abstol) {
8831567311SBarry Smith     mg->ttol = abstol;
894b9ad928SBarry Smith   } else {
9031567311SBarry Smith     mg->ttol = 0.0;
914b9ad928SBarry Smith   }
924b9ad928SBarry Smith 
934d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
944d0a8057SBarry Smith      stop prematurely do to small residual */
954d0a8057SBarry Smith   for (i=1; i<levels; i++) {
96f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
97f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
98f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
994b9ad928SBarry Smith     }
1004d0a8057SBarry Smith   }
1014d0a8057SBarry Smith 
1024d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1034d0a8057SBarry Smith   for (i=0; i<its; i++) {
104f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1054d0a8057SBarry Smith     if (*reason) break;
1064d0a8057SBarry Smith   }
1074d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1084d0a8057SBarry Smith   *outits = i;
1094b9ad928SBarry Smith   PetscFunctionReturn(0);
1104b9ad928SBarry Smith }
1114b9ad928SBarry Smith 
1124b9ad928SBarry Smith #undef __FUNCT__
1139dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
1144b9ad928SBarry Smith /*@C
11597177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1164b9ad928SBarry Smith    Must be called before any other MG routine.
1174b9ad928SBarry Smith 
118ad4df100SBarry Smith    Logically Collective on PC
1194b9ad928SBarry Smith 
1204b9ad928SBarry Smith    Input Parameters:
1214b9ad928SBarry Smith +  pc - the preconditioner context
1224b9ad928SBarry Smith .  levels - the number of levels
1234b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1244b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
1254b9ad928SBarry Smith 
1264b9ad928SBarry Smith    Level: intermediate
1274b9ad928SBarry Smith 
1284b9ad928SBarry Smith    Notes:
1294b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1304b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1314b9ad928SBarry Smith 
1324b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1334b9ad928SBarry Smith 
13497177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1354b9ad928SBarry Smith @*/
1367087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1374b9ad928SBarry Smith {
138dfbe8321SBarry Smith   PetscErrorCode ierr;
139f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
140f3fbd535SBarry Smith   MPI_Comm       comm = ((PetscObject)pc)->comm;
141f3fbd535SBarry Smith   PC_MG_Levels   **mglevels;
142f3fbd535SBarry Smith   PetscInt       i;
143f3fbd535SBarry Smith   PetscMPIInt    size;
144f3fbd535SBarry Smith   const char     *prefix;
145f3fbd535SBarry Smith   PC             ipc;
1464b9ad928SBarry Smith 
1474b9ad928SBarry Smith   PetscFunctionBegin;
1480700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
149*548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
150*548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
151e7e72b3dSBarry 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()");
152e32f2f54SBarry Smith   if (mg->levels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc, this array should not yet exist");
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__
215a06653b4SBarry Smith #define __FUNCT__ "PCReset_MG"
216a06653b4SBarry Smith PetscErrorCode PCReset_MG(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) {
236a06653b4SBarry Smith 	ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
237f3fbd535SBarry Smith       }
238a06653b4SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
239f3fbd535SBarry Smith     }
240f3fbd535SBarry Smith   }
241c07bf074SBarry Smith   PetscFunctionReturn(0);
242c07bf074SBarry Smith }
243c07bf074SBarry Smith 
244c07bf074SBarry Smith #undef __FUNCT__
245c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
246c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
247c07bf074SBarry Smith {
248c07bf074SBarry Smith   PetscErrorCode ierr;
249a06653b4SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
250a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
251a06653b4SBarry Smith   PetscInt       i,n;
252c07bf074SBarry Smith 
253c07bf074SBarry Smith   PetscFunctionBegin;
254a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
255a06653b4SBarry Smith   if (mglevels) {
256a06653b4SBarry Smith     n = mglevels[0]->levels;
257a06653b4SBarry Smith     for (i=0; i<n; i++) {
258a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
259a06653b4SBarry Smith 	ierr = KSPDestroy(mglevels[i]->smoothd);CHKERRQ(ierr);
260a06653b4SBarry Smith       }
261a06653b4SBarry Smith       ierr = KSPDestroy(mglevels[i]->smoothu);CHKERRQ(ierr);
262a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
263a06653b4SBarry Smith     }
264a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
265a06653b4SBarry Smith   }
266c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
267f3fbd535SBarry Smith   PetscFunctionReturn(0);
268f3fbd535SBarry Smith }
269f3fbd535SBarry Smith 
270f3fbd535SBarry Smith 
271f3fbd535SBarry Smith 
27209573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
27309573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
27409573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
275f3fbd535SBarry Smith 
276f3fbd535SBarry Smith /*
277f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
278f3fbd535SBarry Smith              or full cycle of multigrid.
279f3fbd535SBarry Smith 
280f3fbd535SBarry Smith   Note:
281f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
282f3fbd535SBarry Smith */
283f3fbd535SBarry Smith #undef __FUNCT__
284f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
285f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
286f3fbd535SBarry Smith {
287f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
288f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
289f3fbd535SBarry Smith   PetscErrorCode ierr;
290f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
291f3fbd535SBarry Smith 
292f3fbd535SBarry Smith   PetscFunctionBegin;
293e1d8e5deSBarry Smith 
294e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
295e1d8e5deSBarry Smith   for (i=0; i<levels-1; i++) {
296e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
297e1d8e5deSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
298e1d8e5deSBarry Smith     }
299e1d8e5deSBarry Smith   }
300e1d8e5deSBarry Smith 
301f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
302f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
30331567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
304f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
30531567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
306f3fbd535SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
307f3fbd535SBarry Smith     }
308f3fbd535SBarry Smith   }
30931567311SBarry Smith   else if (mg->am == PC_MG_ADDITIVE) {
31031567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
311f3fbd535SBarry Smith   }
31231567311SBarry Smith   else if (mg->am == PC_MG_KASKADE) {
31331567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
314f3fbd535SBarry Smith   }
315f3fbd535SBarry Smith   else {
316f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
317f3fbd535SBarry Smith   }
318f3fbd535SBarry Smith   PetscFunctionReturn(0);
319f3fbd535SBarry Smith }
320f3fbd535SBarry Smith 
321f3fbd535SBarry Smith 
322f3fbd535SBarry Smith #undef __FUNCT__
323f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
324f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc)
325f3fbd535SBarry Smith {
326f3fbd535SBarry Smith   PetscErrorCode ierr;
327f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
328ace3abfcSBarry Smith   PetscBool      flg;
329f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
330f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
331f3fbd535SBarry Smith   PCMGType       mgtype;
332f3fbd535SBarry Smith   PCMGCycleType  mgctype;
333f3fbd535SBarry Smith 
334f3fbd535SBarry Smith   PetscFunctionBegin;
335f3fbd535SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
33618aabeadSBarry Smith     if (!mglevels) {
337f3fbd535SBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
338f3fbd535SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
339f3fbd535SBarry Smith       mglevels = mg->levels;
340f3fbd535SBarry Smith     }
341f3fbd535SBarry Smith     mgctype = (PCMGCycleType) mglevels[0]->cycles;
342f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
343f3fbd535SBarry Smith     if (flg) {
344f3fbd535SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
345f3fbd535SBarry Smith     };
346f3fbd535SBarry Smith     flg  = PETSC_FALSE;
347acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
348f3fbd535SBarry Smith     if (flg) {
349f3fbd535SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
350f3fbd535SBarry Smith     }
351f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
352f3fbd535SBarry Smith     if (flg) {
353f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
354f3fbd535SBarry Smith     }
355f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
356f3fbd535SBarry Smith     if (flg) {
357f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
358f3fbd535SBarry Smith     }
35931567311SBarry Smith     mgtype = mg->am;
360f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
361f3fbd535SBarry Smith     if (flg) {
362f3fbd535SBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
363f3fbd535SBarry Smith     }
36431567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
36531567311SBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
366f3fbd535SBarry Smith       if (flg) {
367f3fbd535SBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
368f3fbd535SBarry Smith       }
369f3fbd535SBarry Smith     }
370f3fbd535SBarry Smith     flg  = PETSC_FALSE;
371acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
372f3fbd535SBarry Smith     if (flg) {
373f3fbd535SBarry Smith       PetscInt i;
374f3fbd535SBarry Smith       char     eventname[128];
37563e6d426SJed Brown       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
376f3fbd535SBarry Smith       levels = mglevels[0]->levels;
377f3fbd535SBarry Smith       for (i=0; i<levels; i++) {
378f3fbd535SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
37963e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
380f3fbd535SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
38163e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
382f3fbd535SBarry Smith         if (i) {
383f3fbd535SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
38463e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
385f3fbd535SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
38663e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
387f3fbd535SBarry Smith         }
388f3fbd535SBarry Smith       }
389f3fbd535SBarry Smith     }
390f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
391f3fbd535SBarry Smith   PetscFunctionReturn(0);
392f3fbd535SBarry Smith }
393f3fbd535SBarry Smith 
394f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
395f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
396f3fbd535SBarry Smith 
397f3fbd535SBarry Smith #undef __FUNCT__
398f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
399f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
400f3fbd535SBarry Smith {
401f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
402f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
403f3fbd535SBarry Smith   PetscErrorCode ierr;
404f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
405ace3abfcSBarry Smith   PetscBool      iascii;
406f3fbd535SBarry Smith 
407f3fbd535SBarry Smith   PetscFunctionBegin;
4082692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
409f3fbd535SBarry Smith   if (iascii) {
41031567311SBarry 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);
41131567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
41231567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
413f3fbd535SBarry Smith     }
414218a07d4SBarry Smith     if (mg->galerkin) {
415f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4164f66f45eSBarry Smith     } else {
4174f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
418f3fbd535SBarry Smith     }
419f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
420f3fbd535SBarry Smith       if (!i) {
4212d19a89fSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
422f3fbd535SBarry Smith       } else {
42331567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
424f3fbd535SBarry Smith       }
425f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
426f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
427f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
428f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
429f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
430f3fbd535SBarry Smith       } else if (i){
43131567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr);
432f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
433f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
434f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
435f3fbd535SBarry Smith       }
436f3fbd535SBarry Smith     }
437f3fbd535SBarry Smith   } else {
43865e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
439f3fbd535SBarry Smith   }
440f3fbd535SBarry Smith   PetscFunctionReturn(0);
441f3fbd535SBarry Smith }
442f3fbd535SBarry Smith 
443f3fbd535SBarry Smith /*
444f3fbd535SBarry Smith     Calls setup for the KSP on each level
445f3fbd535SBarry Smith */
446f3fbd535SBarry Smith #undef __FUNCT__
447f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
448f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
449f3fbd535SBarry Smith {
450f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
451f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
452f3fbd535SBarry Smith   PetscErrorCode          ierr;
453f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
454f3fbd535SBarry Smith   PC                      cpc,mpc;
455ace3abfcSBarry Smith   PetscBool               preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset;
456f3fbd535SBarry Smith   PetscViewerASCIIMonitor ascii;
457f3fbd535SBarry Smith   PetscViewer             viewer = PETSC_NULL;
458f3fbd535SBarry Smith   MPI_Comm                comm;
459f3fbd535SBarry Smith   Mat                     dA,dB;
460f3fbd535SBarry Smith   MatStructure            uflag;
461f3fbd535SBarry Smith   Vec                     tvec;
462218a07d4SBarry Smith   DM                      *dms;
463f3fbd535SBarry Smith 
464f3fbd535SBarry Smith   PetscFunctionBegin;
465f3fbd535SBarry Smith 
466f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
467f3fbd535SBarry Smith   /* so use those from global PC */
468f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
469f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
470f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
471f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr);
4721338a6b9SJed Brown   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) {
473f3fbd535SBarry 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);
474f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
475f3fbd535SBarry Smith   }
476f3fbd535SBarry Smith 
4772d2b81a6SBarry Smith   if (pc->dm && !pc->setupcalled) {
4782d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
4792e499ae9SBarry Smith     Mat p;
480218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
481218a07d4SBarry Smith     dms[n-1] = pc->dm;
482218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
483218a07d4SBarry Smith       ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr);
48411629dbeSBarry Smith       ierr = DMSetFunction(dms[i],0);
48511629dbeSBarry Smith       ierr = DMSetInitialGuess(dms[i],0);
48624c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
4872d2b81a6SBarry Smith 	ierr = DMGetInterpolation(dms[i],dms[i+1],&p,PETSC_NULL);CHKERRQ(ierr);
4882d2b81a6SBarry Smith 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
4892d2b81a6SBarry Smith         ierr = MatDestroy(p);CHKERRQ(ierr);
490218a07d4SBarry Smith       }
49124c3aa18SBarry Smith     }
4922d2b81a6SBarry Smith 
4932d2b81a6SBarry Smith     if (!mg->galerkin) {
494d42688cbSBarry Smith       /* each coarse level gets its DM; finest level does not get DM because it shared the outer PC operators */
4952d2b81a6SBarry Smith       for (i=n-2; i>-1; i--) {
4962e499ae9SBarry Smith         ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
4972d2b81a6SBarry Smith       }
4982d2b81a6SBarry Smith     }
4992d2b81a6SBarry Smith 
500218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
501218a07d4SBarry Smith       ierr = DMDestroy(dms[i]);CHKERRQ(ierr);
502218a07d4SBarry Smith     }
5032d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
504218a07d4SBarry Smith   }
505218a07d4SBarry Smith 
506218a07d4SBarry Smith   if (mg->galerkin) {
507f3fbd535SBarry Smith     Mat B;
508218a07d4SBarry Smith     mg->galerkinused = PETSC_TRUE;
509f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
510f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
511f3fbd535SBarry Smith     if (!pc->setupcalled) {
512f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
513f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
514f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
515f3fbd535SBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
516f3fbd535SBarry Smith         dB   = B;
517f3fbd535SBarry Smith       }
518cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
519f3fbd535SBarry Smith     } else {
520f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
521f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
522f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
523f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
524f3fbd535SBarry Smith         dB   = B;
525f3fbd535SBarry Smith       }
526f3fbd535SBarry Smith     }
527f3fbd535SBarry Smith   }
528f3fbd535SBarry Smith 
529f3fbd535SBarry Smith   if (!pc->setupcalled) {
530acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr);
531f3fbd535SBarry Smith 
532f3fbd535SBarry Smith     for (i=0; i<n; i++) {
533f3fbd535SBarry Smith       if (monitor) {
534f3fbd535SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr);
535f3fbd535SBarry Smith         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
536f3fbd535SBarry Smith         ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
537f3fbd535SBarry Smith       }
538f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
539f3fbd535SBarry Smith     }
540f3fbd535SBarry Smith     for (i=1; i<n; i++) {
541f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
542f3fbd535SBarry Smith         if (monitor) {
543f3fbd535SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr);
544f3fbd535SBarry Smith           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
545f3fbd535SBarry Smith           ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
546f3fbd535SBarry Smith         }
547f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
548f3fbd535SBarry Smith       }
549f3fbd535SBarry Smith     }
550f3fbd535SBarry Smith     for (i=1; i<n; i++) {
551f3fbd535SBarry Smith       if (mglevels[i]->restrct && !mglevels[i]->interpolate) {
552f3fbd535SBarry Smith         ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr);
553f3fbd535SBarry Smith       }
554f3fbd535SBarry Smith       if (!mglevels[i]->restrct && mglevels[i]->interpolate) {
555f3fbd535SBarry Smith         ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr);
556f3fbd535SBarry Smith       }
557f3fbd535SBarry Smith #if defined(PETSC_USE_DEBUG)
558f3fbd535SBarry Smith       if (!mglevels[i]->restrct || !mglevels[i]->interpolate) {
55965e19b50SBarry Smith         SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
560f3fbd535SBarry Smith       }
561f3fbd535SBarry Smith #endif
562f3fbd535SBarry Smith     }
563f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
564f3fbd535SBarry Smith       if (!mglevels[i]->b) {
565f3fbd535SBarry Smith         Vec *vec;
566f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
567f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
568f3fbd535SBarry Smith         ierr = VecDestroy(*vec);CHKERRQ(ierr);
569f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
570f3fbd535SBarry Smith       }
571f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
572f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
573f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
574f3fbd535SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
575f3fbd535SBarry Smith       }
576f3fbd535SBarry Smith       if (!mglevels[i]->x) {
577f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
578f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
579f3fbd535SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
580f3fbd535SBarry Smith       }
581f3fbd535SBarry Smith     }
582f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
583f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
584f3fbd535SBarry Smith       Vec *vec;
585f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
586f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
587f3fbd535SBarry Smith       ierr = VecDestroy(*vec);CHKERRQ(ierr);
588f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
589f3fbd535SBarry Smith     }
590f3fbd535SBarry Smith   }
591f3fbd535SBarry Smith 
592f3fbd535SBarry Smith 
593f3fbd535SBarry Smith   for (i=1; i<n; i++) {
594f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
595f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
596f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
597f3fbd535SBarry Smith     }
59863e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
599f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
60063e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
601d42688cbSBarry Smith     if (!mglevels[i]->residual) {
602d42688cbSBarry Smith       Mat mat;
603d42688cbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
604d42688cbSBarry Smith       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
605d42688cbSBarry Smith     }
606f3fbd535SBarry Smith   }
607f3fbd535SBarry Smith   for (i=1; i<n; i++) {
608f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
609f3fbd535SBarry Smith       Mat          downmat,downpmat;
610f3fbd535SBarry Smith       MatStructure matflag;
611ace3abfcSBarry Smith       PetscBool    opsset;
612f3fbd535SBarry Smith 
613f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
614f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
615f3fbd535SBarry Smith       if (!opsset) {
616f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
617f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
618f3fbd535SBarry Smith       }
619f3fbd535SBarry Smith 
620f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
62163e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
622f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
62363e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
624f3fbd535SBarry Smith     }
625f3fbd535SBarry Smith   }
626f3fbd535SBarry Smith 
627f3fbd535SBarry Smith   /*
628f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
629f3fbd535SBarry Smith   */
630f3fbd535SBarry Smith   ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
631f3fbd535SBarry Smith   if (preonly) {
632f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
633f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
634f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
635f3fbd535SBarry Smith     if (!lu && !redundant && !cholesky) {
636f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
637f3fbd535SBarry Smith     }
638f3fbd535SBarry Smith   }
639f3fbd535SBarry Smith 
640f3fbd535SBarry Smith   if (!pc->setupcalled) {
641f3fbd535SBarry Smith     if (monitor) {
642f3fbd535SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr);
643f3fbd535SBarry Smith       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
644f3fbd535SBarry Smith       ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
645f3fbd535SBarry Smith     }
646f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
647f3fbd535SBarry Smith   }
648f3fbd535SBarry Smith 
64963e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
650f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
65163e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
652f3fbd535SBarry Smith 
653f3fbd535SBarry Smith   /*
654f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
655e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
656f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
657f3fbd535SBarry Smith 
658f3fbd535SBarry Smith    Only support one or the other at the same time.
659f3fbd535SBarry Smith   */
660f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
661acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
662f3fbd535SBarry Smith   if (dump) {
663f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
664f3fbd535SBarry Smith   }
665f3fbd535SBarry Smith   dump = PETSC_FALSE;
666f3fbd535SBarry Smith #endif
667acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
668f3fbd535SBarry Smith   if (dump) {
669f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
670f3fbd535SBarry Smith   }
671f3fbd535SBarry Smith 
672f3fbd535SBarry Smith   if (viewer) {
673f3fbd535SBarry Smith     for (i=1; i<n; i++) {
674f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
675f3fbd535SBarry Smith     }
676f3fbd535SBarry Smith     for (i=0; i<n; i++) {
677f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
678f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
679f3fbd535SBarry Smith     }
680f3fbd535SBarry Smith   }
681f3fbd535SBarry Smith   PetscFunctionReturn(0);
682f3fbd535SBarry Smith }
683f3fbd535SBarry Smith 
684f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
685f3fbd535SBarry Smith 
686f3fbd535SBarry Smith #undef __FUNCT__
6879dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
6884b9ad928SBarry Smith /*@
68997177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
6904b9ad928SBarry Smith 
6914b9ad928SBarry Smith    Not Collective
6924b9ad928SBarry Smith 
6934b9ad928SBarry Smith    Input Parameter:
6944b9ad928SBarry Smith .  pc - the preconditioner context
6954b9ad928SBarry Smith 
6964b9ad928SBarry Smith    Output parameter:
6974b9ad928SBarry Smith .  levels - the number of levels
6984b9ad928SBarry Smith 
6994b9ad928SBarry Smith    Level: advanced
7004b9ad928SBarry Smith 
7014b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
7024b9ad928SBarry Smith 
70397177400SBarry Smith .seealso: PCMGSetLevels()
7044b9ad928SBarry Smith @*/
7057087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
7064b9ad928SBarry Smith {
707f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
7084b9ad928SBarry Smith 
7094b9ad928SBarry Smith   PetscFunctionBegin;
7100700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
7114482741eSBarry Smith   PetscValidIntPointer(levels,2);
712f3fbd535SBarry Smith   *levels = mg->nlevels;
7134b9ad928SBarry Smith   PetscFunctionReturn(0);
7144b9ad928SBarry Smith }
7154b9ad928SBarry Smith 
7164b9ad928SBarry Smith #undef __FUNCT__
7179dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
7184b9ad928SBarry Smith /*@
71997177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
7204b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
7214b9ad928SBarry Smith 
722ad4df100SBarry Smith    Logically Collective on PC
7234b9ad928SBarry Smith 
7244b9ad928SBarry Smith    Input Parameters:
7254b9ad928SBarry Smith +  pc - the preconditioner context
7269dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
7279dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
7284b9ad928SBarry Smith 
7294b9ad928SBarry Smith    Options Database Key:
7304b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
7314b9ad928SBarry Smith    additive, full, kaskade
7324b9ad928SBarry Smith 
7334b9ad928SBarry Smith    Level: advanced
7344b9ad928SBarry Smith 
7354b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
7364b9ad928SBarry Smith 
73797177400SBarry Smith .seealso: PCMGSetLevels()
7384b9ad928SBarry Smith @*/
7397087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
7404b9ad928SBarry Smith {
741f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
7424b9ad928SBarry Smith 
7434b9ad928SBarry Smith   PetscFunctionBegin;
7440700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
745c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
74631567311SBarry Smith   mg->am = form;
7479dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
7484b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
7494b9ad928SBarry Smith   PetscFunctionReturn(0);
7504b9ad928SBarry Smith }
7514b9ad928SBarry Smith 
7524b9ad928SBarry Smith #undef __FUNCT__
7530d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
7544b9ad928SBarry Smith /*@
7550d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
7564b9ad928SBarry Smith    complicated cycling.
7574b9ad928SBarry Smith 
758ad4df100SBarry Smith    Logically Collective on PC
7594b9ad928SBarry Smith 
7604b9ad928SBarry Smith    Input Parameters:
761c2be2410SBarry Smith +  pc - the multigrid context
7620d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
7634b9ad928SBarry Smith 
7644b9ad928SBarry Smith    Options Database Key:
7650d353602SBarry Smith $  -pc_mg_cycle_type v or w
7664b9ad928SBarry Smith 
7674b9ad928SBarry Smith    Level: advanced
7684b9ad928SBarry Smith 
7694b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7704b9ad928SBarry Smith 
7710d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
7724b9ad928SBarry Smith @*/
7737087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
7744b9ad928SBarry Smith {
775f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
776f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
77779416396SBarry Smith   PetscInt     i,levels;
7784b9ad928SBarry Smith 
7794b9ad928SBarry Smith   PetscFunctionBegin;
7800700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
781e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
782c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
783f3fbd535SBarry Smith   levels = mglevels[0]->levels;
7844b9ad928SBarry Smith 
7854b9ad928SBarry Smith   for (i=0; i<levels; i++) {
786f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
7874b9ad928SBarry Smith   }
7884b9ad928SBarry Smith   PetscFunctionReturn(0);
7894b9ad928SBarry Smith }
7904b9ad928SBarry Smith 
7914b9ad928SBarry Smith #undef __FUNCT__
7928cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
7938cc2d5dfSBarry Smith /*@
7948cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
7958cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
7968cc2d5dfSBarry Smith 
797ad4df100SBarry Smith    Logically Collective on PC
7988cc2d5dfSBarry Smith 
7998cc2d5dfSBarry Smith    Input Parameters:
8008cc2d5dfSBarry Smith +  pc - the multigrid context
8018cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
8028cc2d5dfSBarry Smith 
8038cc2d5dfSBarry Smith    Options Database Key:
8048cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
8058cc2d5dfSBarry Smith 
8068cc2d5dfSBarry Smith    Level: advanced
8078cc2d5dfSBarry Smith 
8088cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
8098cc2d5dfSBarry Smith 
8108cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8118cc2d5dfSBarry Smith 
8128cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
8138cc2d5dfSBarry Smith @*/
8147087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
8158cc2d5dfSBarry Smith {
816f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
817f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
8188cc2d5dfSBarry Smith   PetscInt     i,levels;
8198cc2d5dfSBarry Smith 
8208cc2d5dfSBarry Smith   PetscFunctionBegin;
8210700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
822e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
823c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
824f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8258cc2d5dfSBarry Smith 
8268cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
82731567311SBarry Smith     mg->cyclesperpcapply  = n;
8288cc2d5dfSBarry Smith   }
8298cc2d5dfSBarry Smith   PetscFunctionReturn(0);
8308cc2d5dfSBarry Smith }
8318cc2d5dfSBarry Smith 
8328cc2d5dfSBarry Smith #undef __FUNCT__
8339dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
834c2be2410SBarry Smith /*@
83597177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
836c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
837c2be2410SBarry Smith 
838ad4df100SBarry Smith    Logically Collective on PC
839c2be2410SBarry Smith 
840c2be2410SBarry Smith    Input Parameters:
8413fc8bf9cSBarry Smith .  pc - the multigrid context
842c2be2410SBarry Smith 
843c2be2410SBarry Smith    Options Database Key:
844c2be2410SBarry Smith $  -pc_mg_galerkin
845c2be2410SBarry Smith 
846c2be2410SBarry Smith    Level: intermediate
847c2be2410SBarry Smith 
848c2be2410SBarry Smith .keywords: MG, set, Galerkin
849c2be2410SBarry Smith 
8503fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
8513fc8bf9cSBarry Smith 
852c2be2410SBarry Smith @*/
8537087cfbeSBarry Smith PetscErrorCode  PCMGSetGalerkin(PC pc)
854c2be2410SBarry Smith {
855f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
856c2be2410SBarry Smith 
857c2be2410SBarry Smith   PetscFunctionBegin;
8580700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
859218a07d4SBarry Smith   mg->galerkin = PETSC_TRUE;
860c2be2410SBarry Smith   PetscFunctionReturn(0);
861c2be2410SBarry Smith }
862c2be2410SBarry Smith 
863c2be2410SBarry Smith #undef __FUNCT__
8643fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
8653fc8bf9cSBarry Smith /*@
8663fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
8673fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
8683fc8bf9cSBarry Smith 
8693fc8bf9cSBarry Smith    Not Collective
8703fc8bf9cSBarry Smith 
8713fc8bf9cSBarry Smith    Input Parameter:
8723fc8bf9cSBarry Smith .  pc - the multigrid context
8733fc8bf9cSBarry Smith 
8743fc8bf9cSBarry Smith    Output Parameter:
8753fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
8763fc8bf9cSBarry Smith 
8773fc8bf9cSBarry Smith    Options Database Key:
8783fc8bf9cSBarry Smith $  -pc_mg_galerkin
8793fc8bf9cSBarry Smith 
8803fc8bf9cSBarry Smith    Level: intermediate
8813fc8bf9cSBarry Smith 
8823fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
8833fc8bf9cSBarry Smith 
8843fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
8853fc8bf9cSBarry Smith 
8863fc8bf9cSBarry Smith @*/
8877087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
8883fc8bf9cSBarry Smith {
889f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
8903fc8bf9cSBarry Smith 
8913fc8bf9cSBarry Smith   PetscFunctionBegin;
8920700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
893218a07d4SBarry Smith   *galerkin = mg->galerkin;
8943fc8bf9cSBarry Smith   PetscFunctionReturn(0);
8953fc8bf9cSBarry Smith }
8963fc8bf9cSBarry Smith 
8973fc8bf9cSBarry Smith #undef __FUNCT__
8989dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
8994b9ad928SBarry Smith /*@
90097177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
90197177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
9024b9ad928SBarry Smith    pre-smoothing steps on different levels.
9034b9ad928SBarry Smith 
904ad4df100SBarry Smith    Logically Collective on PC
9054b9ad928SBarry Smith 
9064b9ad928SBarry Smith    Input Parameters:
9074b9ad928SBarry Smith +  mg - the multigrid context
9084b9ad928SBarry Smith -  n - the number of smoothing steps
9094b9ad928SBarry Smith 
9104b9ad928SBarry Smith    Options Database Key:
9114b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
9124b9ad928SBarry Smith 
9134b9ad928SBarry Smith    Level: advanced
9144b9ad928SBarry Smith 
9154b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
9164b9ad928SBarry Smith 
91797177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
9184b9ad928SBarry Smith @*/
9197087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
9204b9ad928SBarry Smith {
921f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
922f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9236849ba73SBarry Smith   PetscErrorCode ierr;
92479416396SBarry Smith   PetscInt       i,levels;
9254b9ad928SBarry Smith 
9264b9ad928SBarry Smith   PetscFunctionBegin;
9270700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
928e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
929c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
930f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9314b9ad928SBarry Smith 
932b05257ddSBarry Smith   for (i=1; i<levels; i++) {
9334b9ad928SBarry Smith     /* make sure smoother up and down are different */
93497177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
935f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
93631567311SBarry Smith     mg->default_smoothd = n;
9374b9ad928SBarry Smith   }
9384b9ad928SBarry Smith   PetscFunctionReturn(0);
9394b9ad928SBarry Smith }
9404b9ad928SBarry Smith 
9414b9ad928SBarry Smith #undef __FUNCT__
9429dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
9434b9ad928SBarry Smith /*@
94497177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
94597177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
9464b9ad928SBarry Smith    post-smoothing steps on different levels.
9474b9ad928SBarry Smith 
948ad4df100SBarry Smith    Logically Collective on PC
9494b9ad928SBarry Smith 
9504b9ad928SBarry Smith    Input Parameters:
9514b9ad928SBarry Smith +  mg - the multigrid context
9524b9ad928SBarry Smith -  n - the number of smoothing steps
9534b9ad928SBarry Smith 
9544b9ad928SBarry Smith    Options Database Key:
9554b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
9564b9ad928SBarry Smith 
9574b9ad928SBarry Smith    Level: advanced
9584b9ad928SBarry Smith 
9594b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
960a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
9614b9ad928SBarry Smith 
9624b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
9634b9ad928SBarry Smith 
96497177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
9654b9ad928SBarry Smith @*/
9667087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
9674b9ad928SBarry Smith {
968f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
969f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9706849ba73SBarry Smith   PetscErrorCode ierr;
97179416396SBarry Smith   PetscInt       i,levels;
9724b9ad928SBarry Smith 
9734b9ad928SBarry Smith   PetscFunctionBegin;
9740700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
975e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
976c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
977f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9784b9ad928SBarry Smith 
9794b9ad928SBarry Smith   for (i=1; i<levels; i++) {
9804b9ad928SBarry Smith     /* make sure smoother up and down are different */
98197177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
982f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
98331567311SBarry Smith     mg->default_smoothu = n;
9844b9ad928SBarry Smith   }
9854b9ad928SBarry Smith   PetscFunctionReturn(0);
9864b9ad928SBarry Smith }
9874b9ad928SBarry Smith 
9884b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
9894b9ad928SBarry Smith 
9903b09bd56SBarry Smith /*MC
991ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
9923b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
9933b09bd56SBarry Smith 
9943b09bd56SBarry Smith    Options Database Keys:
9953b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
9960d353602SBarry Smith .  -pc_mg_cycles v or w
99779416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
9983b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
9993b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
10003b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
10013b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
100268eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
10033b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1004e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
10053b09bd56SBarry Smith 
100624c3aa18SBarry 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
10073b09bd56SBarry Smith 
10083b09bd56SBarry Smith    Level: intermediate
10093b09bd56SBarry Smith 
10108f87f92bSBarry Smith    Concepts: multigrid/multilevel
10113b09bd56SBarry Smith 
101224c3aa18SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC
10130d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
101497177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
101597177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
10160d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
10173b09bd56SBarry Smith M*/
10183b09bd56SBarry Smith 
10194b9ad928SBarry Smith EXTERN_C_BEGIN
10204b9ad928SBarry Smith #undef __FUNCT__
10214b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
10227087cfbeSBarry Smith PetscErrorCode  PCCreate_MG(PC pc)
10234b9ad928SBarry Smith {
1024f3fbd535SBarry Smith   PC_MG          *mg;
1025f3fbd535SBarry Smith   PetscErrorCode ierr;
1026f3fbd535SBarry Smith 
10274b9ad928SBarry Smith   PetscFunctionBegin;
1028f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1029f3fbd535SBarry Smith   pc->data    = (void*)mg;
1030f3fbd535SBarry Smith   mg->nlevels = -1;
1031f3fbd535SBarry Smith 
10324b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
10334b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1034a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
10354b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
10364b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
10374b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
10384b9ad928SBarry Smith   PetscFunctionReturn(0);
10394b9ad928SBarry Smith }
10404b9ad928SBarry Smith EXTERN_C_END
1041