xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 3c0c59f39ee3bd561a1b4dcbd0e58366bfa60aba)
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);
149548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
150548767e0SJed 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++) {
2276bf464f9SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
2286bf464f9SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
2296bf464f9SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
2306bf464f9SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
2316bf464f9SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
232f3fbd535SBarry Smith     }
233f3fbd535SBarry Smith 
234f3fbd535SBarry Smith     for (i=0; i<n; i++) {
235298cc208SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
236f3fbd535SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
237a06653b4SBarry Smith 	ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
238f3fbd535SBarry Smith       }
239a06653b4SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
240f3fbd535SBarry Smith     }
241f3fbd535SBarry Smith   }
242c07bf074SBarry Smith   PetscFunctionReturn(0);
243c07bf074SBarry Smith }
244c07bf074SBarry Smith 
245c07bf074SBarry Smith #undef __FUNCT__
246c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
247c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
248c07bf074SBarry Smith {
249c07bf074SBarry Smith   PetscErrorCode ierr;
250a06653b4SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
251a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
252a06653b4SBarry Smith   PetscInt       i,n;
253c07bf074SBarry Smith 
254c07bf074SBarry Smith   PetscFunctionBegin;
255a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
256a06653b4SBarry Smith   if (mglevels) {
257a06653b4SBarry Smith     n = mglevels[0]->levels;
258a06653b4SBarry Smith     for (i=0; i<n; i++) {
259a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2606bf464f9SBarry Smith 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
261a06653b4SBarry Smith       }
2626bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
263a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
264a06653b4SBarry Smith     }
265a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
266a06653b4SBarry Smith   }
267c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
268f3fbd535SBarry Smith   PetscFunctionReturn(0);
269f3fbd535SBarry Smith }
270f3fbd535SBarry Smith 
271f3fbd535SBarry Smith 
272f3fbd535SBarry Smith 
27309573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
27409573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
27509573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
276f3fbd535SBarry Smith 
277f3fbd535SBarry Smith /*
278f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
279f3fbd535SBarry Smith              or full cycle of multigrid.
280f3fbd535SBarry Smith 
281f3fbd535SBarry Smith   Note:
282f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
283f3fbd535SBarry Smith */
284f3fbd535SBarry Smith #undef __FUNCT__
285f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
286f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
287f3fbd535SBarry Smith {
288f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
289f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
290f3fbd535SBarry Smith   PetscErrorCode ierr;
291f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
292f3fbd535SBarry Smith 
293f3fbd535SBarry Smith   PetscFunctionBegin;
294e1d8e5deSBarry Smith 
295e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
296298cc208SBarry Smith   for (i=0; i<levels; i++) {
297e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
298e1d8e5deSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
299298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
300e1d8e5deSBarry Smith     }
301e1d8e5deSBarry Smith   }
302e1d8e5deSBarry Smith 
303f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
304f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
30531567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
306f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
30731567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
308f3fbd535SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
309f3fbd535SBarry Smith     }
310f3fbd535SBarry Smith   }
31131567311SBarry Smith   else if (mg->am == PC_MG_ADDITIVE) {
31231567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
313f3fbd535SBarry Smith   }
31431567311SBarry Smith   else if (mg->am == PC_MG_KASKADE) {
31531567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
316f3fbd535SBarry Smith   }
317f3fbd535SBarry Smith   else {
318f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
319f3fbd535SBarry Smith   }
320f3fbd535SBarry Smith   PetscFunctionReturn(0);
321f3fbd535SBarry Smith }
322f3fbd535SBarry Smith 
323f3fbd535SBarry Smith 
324f3fbd535SBarry Smith #undef __FUNCT__
325f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
326f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc)
327f3fbd535SBarry Smith {
328f3fbd535SBarry Smith   PetscErrorCode ierr;
329f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
330ace3abfcSBarry Smith   PetscBool      flg;
331f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
332f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
333f3fbd535SBarry Smith   PCMGType       mgtype;
334f3fbd535SBarry Smith   PCMGCycleType  mgctype;
335f3fbd535SBarry Smith 
336f3fbd535SBarry Smith   PetscFunctionBegin;
337f3fbd535SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
33818aabeadSBarry Smith     if (!mglevels) {
339f3fbd535SBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
340f3fbd535SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
341f3fbd535SBarry Smith       mglevels = mg->levels;
342f3fbd535SBarry Smith     }
343f3fbd535SBarry Smith     mgctype = (PCMGCycleType) mglevels[0]->cycles;
344f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
345f3fbd535SBarry Smith     if (flg) {
346f3fbd535SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
347f3fbd535SBarry Smith     };
348f3fbd535SBarry Smith     flg  = PETSC_FALSE;
349acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
350f3fbd535SBarry Smith     if (flg) {
351f3fbd535SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
352f3fbd535SBarry Smith     }
353f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
354f3fbd535SBarry Smith     if (flg) {
355f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
356f3fbd535SBarry Smith     }
357f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
358f3fbd535SBarry Smith     if (flg) {
359f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
360f3fbd535SBarry Smith     }
36131567311SBarry Smith     mgtype = mg->am;
362f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
363f3fbd535SBarry Smith     if (flg) {
364f3fbd535SBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
365f3fbd535SBarry Smith     }
36631567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
36731567311SBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
368f3fbd535SBarry Smith       if (flg) {
369f3fbd535SBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
370f3fbd535SBarry Smith       }
371f3fbd535SBarry Smith     }
372f3fbd535SBarry Smith     flg  = PETSC_FALSE;
373acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
374f3fbd535SBarry Smith     if (flg) {
375f3fbd535SBarry Smith       PetscInt i;
376f3fbd535SBarry Smith       char     eventname[128];
37763e6d426SJed Brown       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
378f3fbd535SBarry Smith       levels = mglevels[0]->levels;
379f3fbd535SBarry Smith       for (i=0; i<levels; i++) {
380f3fbd535SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
38163e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
382f3fbd535SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
38363e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
384f3fbd535SBarry Smith         if (i) {
385f3fbd535SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
38663e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
387f3fbd535SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
38863e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
389f3fbd535SBarry Smith         }
390f3fbd535SBarry Smith       }
391f3fbd535SBarry Smith     }
392f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
393f3fbd535SBarry Smith   PetscFunctionReturn(0);
394f3fbd535SBarry Smith }
395f3fbd535SBarry Smith 
396f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
397f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
398f3fbd535SBarry Smith 
399f3fbd535SBarry Smith #undef __FUNCT__
400f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
401f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
402f3fbd535SBarry Smith {
403f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
404f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
405f3fbd535SBarry Smith   PetscErrorCode ierr;
406f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
407ace3abfcSBarry Smith   PetscBool      iascii;
408f3fbd535SBarry Smith 
409f3fbd535SBarry Smith   PetscFunctionBegin;
4102692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
411f3fbd535SBarry Smith   if (iascii) {
41231567311SBarry 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);
41331567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
41431567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
415f3fbd535SBarry Smith     }
416218a07d4SBarry Smith     if (mg->galerkin) {
417f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4184f66f45eSBarry Smith     } else {
4194f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
420f3fbd535SBarry Smith     }
421f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
422f3fbd535SBarry Smith       if (!i) {
4232d19a89fSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
424f3fbd535SBarry Smith       } else {
42531567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
426f3fbd535SBarry Smith       }
427f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
428f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
429f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
430f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
431f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
432f3fbd535SBarry Smith       } else if (i){
43331567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr);
434f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
435f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
436f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
437f3fbd535SBarry Smith       }
438f3fbd535SBarry Smith     }
439f3fbd535SBarry Smith   } else {
44065e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
441f3fbd535SBarry Smith   }
442f3fbd535SBarry Smith   PetscFunctionReturn(0);
443f3fbd535SBarry Smith }
444f3fbd535SBarry Smith 
445f3fbd535SBarry Smith /*
446f3fbd535SBarry Smith     Calls setup for the KSP on each level
447f3fbd535SBarry Smith */
448f3fbd535SBarry Smith #undef __FUNCT__
449f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
450f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
451f3fbd535SBarry Smith {
452f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
453f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
454f3fbd535SBarry Smith   PetscErrorCode          ierr;
455f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
456f3fbd535SBarry Smith   PC                      cpc,mpc;
457ace3abfcSBarry Smith   PetscBool               preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset;
458f3fbd535SBarry Smith   PetscViewerASCIIMonitor ascii;
459f3fbd535SBarry Smith   PetscViewer             viewer = PETSC_NULL;
460f3fbd535SBarry Smith   MPI_Comm                comm;
461f3fbd535SBarry Smith   Mat                     dA,dB;
462f3fbd535SBarry Smith   MatStructure            uflag;
463f3fbd535SBarry Smith   Vec                     tvec;
464218a07d4SBarry Smith   DM                      *dms;
465f3fbd535SBarry Smith 
466f3fbd535SBarry Smith   PetscFunctionBegin;
467f3fbd535SBarry Smith 
468f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
469f3fbd535SBarry Smith   /* so use those from global PC */
470f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
471f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
472f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
473f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr);
4741338a6b9SJed Brown   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) {
475f3fbd535SBarry 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);
476f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
477f3fbd535SBarry Smith   }
478f3fbd535SBarry Smith 
4792d2b81a6SBarry Smith   if (pc->dm && !pc->setupcalled) {
4802d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
4812e499ae9SBarry Smith     Mat p;
482218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
483218a07d4SBarry Smith     dms[n-1] = pc->dm;
484218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
485218a07d4SBarry Smith       ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr);
486*3c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
487*3c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
48811629dbeSBarry Smith       ierr = DMSetFunction(dms[i],0);
48911629dbeSBarry Smith       ierr = DMSetInitialGuess(dms[i],0);
49024c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
4912d2b81a6SBarry Smith 	ierr = DMGetInterpolation(dms[i],dms[i+1],&p,PETSC_NULL);CHKERRQ(ierr);
4922d2b81a6SBarry Smith 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
4936bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
494218a07d4SBarry Smith       }
49524c3aa18SBarry Smith     }
4962d2b81a6SBarry Smith 
497218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
4986bf464f9SBarry Smith       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
499218a07d4SBarry Smith     }
5002d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
501218a07d4SBarry Smith   }
502218a07d4SBarry Smith 
503218a07d4SBarry Smith   if (mg->galerkin) {
504f3fbd535SBarry Smith     Mat B;
505218a07d4SBarry Smith     mg->galerkinused = PETSC_TRUE;
506f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
507f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
508f3fbd535SBarry Smith     if (!pc->setupcalled) {
509f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
510f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
511f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
512f3fbd535SBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
513f3fbd535SBarry Smith         dB   = B;
514f3fbd535SBarry Smith       }
515cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
516f3fbd535SBarry Smith     } else {
517f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
518f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
519f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
520f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
521f3fbd535SBarry Smith         dB   = B;
522f3fbd535SBarry Smith       }
523f3fbd535SBarry Smith     }
524f3fbd535SBarry Smith   }
525f3fbd535SBarry Smith 
526f3fbd535SBarry Smith   if (!pc->setupcalled) {
527acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr);
528f3fbd535SBarry Smith 
529f3fbd535SBarry Smith     for (i=0; i<n; i++) {
530f3fbd535SBarry Smith       if (monitor) {
531f3fbd535SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr);
532f3fbd535SBarry Smith         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
533c2efdce3SBarry Smith         ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
534f3fbd535SBarry Smith       }
535f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
536f3fbd535SBarry Smith     }
537f3fbd535SBarry Smith     for (i=1; i<n; i++) {
538f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
539f3fbd535SBarry Smith         if (monitor) {
540f3fbd535SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr);
541f3fbd535SBarry Smith           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
542c2efdce3SBarry Smith           ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
543f3fbd535SBarry Smith         }
544f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
545f3fbd535SBarry Smith       }
546f3fbd535SBarry Smith     }
547f3fbd535SBarry Smith     for (i=1; i<n; i++) {
548f3fbd535SBarry Smith       if (mglevels[i]->restrct && !mglevels[i]->interpolate) {
549f3fbd535SBarry Smith         ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr);
550f3fbd535SBarry Smith       }
551f3fbd535SBarry Smith       if (!mglevels[i]->restrct && mglevels[i]->interpolate) {
552f3fbd535SBarry Smith         ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr);
553f3fbd535SBarry Smith       }
554f3fbd535SBarry Smith #if defined(PETSC_USE_DEBUG)
555f3fbd535SBarry Smith       if (!mglevels[i]->restrct || !mglevels[i]->interpolate) {
55665e19b50SBarry Smith         SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
557f3fbd535SBarry Smith       }
558f3fbd535SBarry Smith #endif
559f3fbd535SBarry Smith     }
560f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
561f3fbd535SBarry Smith       if (!mglevels[i]->b) {
562f3fbd535SBarry Smith         Vec *vec;
563f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
564f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
5656bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
566f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
567f3fbd535SBarry Smith       }
568f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
569f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
570f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
5716bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
572f3fbd535SBarry Smith       }
573f3fbd535SBarry Smith       if (!mglevels[i]->x) {
574f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
575f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
5766bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
577f3fbd535SBarry Smith       }
578f3fbd535SBarry Smith     }
579f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
580f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
581f3fbd535SBarry Smith       Vec *vec;
582f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
583f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
5846bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
585f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
586f3fbd535SBarry Smith     }
587f3fbd535SBarry Smith   }
588f3fbd535SBarry Smith 
589f3fbd535SBarry Smith 
590f3fbd535SBarry Smith   for (i=1; i<n; i++) {
591f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
592f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
593f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
594f3fbd535SBarry Smith     }
59563e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
596f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
59763e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
598d42688cbSBarry Smith     if (!mglevels[i]->residual) {
599d42688cbSBarry Smith       Mat mat;
600d42688cbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
601d42688cbSBarry Smith       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
602d42688cbSBarry Smith     }
603f3fbd535SBarry Smith   }
604f3fbd535SBarry Smith   for (i=1; i<n; i++) {
605f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
606f3fbd535SBarry Smith       Mat          downmat,downpmat;
607f3fbd535SBarry Smith       MatStructure matflag;
608ace3abfcSBarry Smith       PetscBool    opsset;
609f3fbd535SBarry Smith 
610f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
611f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
612f3fbd535SBarry Smith       if (!opsset) {
613f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
614f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
615f3fbd535SBarry Smith       }
616f3fbd535SBarry Smith 
617f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
61863e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
619f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
62063e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
621f3fbd535SBarry Smith     }
622f3fbd535SBarry Smith   }
623f3fbd535SBarry Smith 
624f3fbd535SBarry Smith   /*
625f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
626f3fbd535SBarry Smith   */
627f3fbd535SBarry Smith   ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
628f3fbd535SBarry Smith   if (preonly) {
629f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
630f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
631f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
632f3fbd535SBarry Smith     if (!lu && !redundant && !cholesky) {
633f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
634f3fbd535SBarry Smith     }
635f3fbd535SBarry Smith   }
636f3fbd535SBarry Smith 
637f3fbd535SBarry Smith   if (!pc->setupcalled) {
638f3fbd535SBarry Smith     if (monitor) {
639f3fbd535SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr);
640f3fbd535SBarry Smith       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
641c2efdce3SBarry Smith       ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
642f3fbd535SBarry Smith     }
643f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
644f3fbd535SBarry Smith   }
645f3fbd535SBarry Smith 
64663e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
647f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
64863e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
649f3fbd535SBarry Smith 
650f3fbd535SBarry Smith   /*
651f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
652e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
653f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
654f3fbd535SBarry Smith 
655f3fbd535SBarry Smith    Only support one or the other at the same time.
656f3fbd535SBarry Smith   */
657f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
658acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
659f3fbd535SBarry Smith   if (dump) {
660f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
661f3fbd535SBarry Smith   }
662f3fbd535SBarry Smith   dump = PETSC_FALSE;
663f3fbd535SBarry Smith #endif
664acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
665f3fbd535SBarry Smith   if (dump) {
666f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
667f3fbd535SBarry Smith   }
668f3fbd535SBarry Smith 
669f3fbd535SBarry Smith   if (viewer) {
670f3fbd535SBarry Smith     for (i=1; i<n; i++) {
671f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
672f3fbd535SBarry Smith     }
673f3fbd535SBarry Smith     for (i=0; i<n; i++) {
674f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
675f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
676f3fbd535SBarry Smith     }
677f3fbd535SBarry Smith   }
678f3fbd535SBarry Smith   PetscFunctionReturn(0);
679f3fbd535SBarry Smith }
680f3fbd535SBarry Smith 
681f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
682f3fbd535SBarry Smith 
683f3fbd535SBarry Smith #undef __FUNCT__
6849dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
6854b9ad928SBarry Smith /*@
68697177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
6874b9ad928SBarry Smith 
6884b9ad928SBarry Smith    Not Collective
6894b9ad928SBarry Smith 
6904b9ad928SBarry Smith    Input Parameter:
6914b9ad928SBarry Smith .  pc - the preconditioner context
6924b9ad928SBarry Smith 
6934b9ad928SBarry Smith    Output parameter:
6944b9ad928SBarry Smith .  levels - the number of levels
6954b9ad928SBarry Smith 
6964b9ad928SBarry Smith    Level: advanced
6974b9ad928SBarry Smith 
6984b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
6994b9ad928SBarry Smith 
70097177400SBarry Smith .seealso: PCMGSetLevels()
7014b9ad928SBarry Smith @*/
7027087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
7034b9ad928SBarry Smith {
704f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
7054b9ad928SBarry Smith 
7064b9ad928SBarry Smith   PetscFunctionBegin;
7070700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
7084482741eSBarry Smith   PetscValidIntPointer(levels,2);
709f3fbd535SBarry Smith   *levels = mg->nlevels;
7104b9ad928SBarry Smith   PetscFunctionReturn(0);
7114b9ad928SBarry Smith }
7124b9ad928SBarry Smith 
7134b9ad928SBarry Smith #undef __FUNCT__
7149dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
7154b9ad928SBarry Smith /*@
71697177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
7174b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
7184b9ad928SBarry Smith 
719ad4df100SBarry Smith    Logically Collective on PC
7204b9ad928SBarry Smith 
7214b9ad928SBarry Smith    Input Parameters:
7224b9ad928SBarry Smith +  pc - the preconditioner context
7239dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
7249dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
7254b9ad928SBarry Smith 
7264b9ad928SBarry Smith    Options Database Key:
7274b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
7284b9ad928SBarry Smith    additive, full, kaskade
7294b9ad928SBarry Smith 
7304b9ad928SBarry Smith    Level: advanced
7314b9ad928SBarry Smith 
7324b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
7334b9ad928SBarry Smith 
73497177400SBarry Smith .seealso: PCMGSetLevels()
7354b9ad928SBarry Smith @*/
7367087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
7374b9ad928SBarry Smith {
738f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
7394b9ad928SBarry Smith 
7404b9ad928SBarry Smith   PetscFunctionBegin;
7410700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
742c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
74331567311SBarry Smith   mg->am = form;
7449dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
7454b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
7464b9ad928SBarry Smith   PetscFunctionReturn(0);
7474b9ad928SBarry Smith }
7484b9ad928SBarry Smith 
7494b9ad928SBarry Smith #undef __FUNCT__
7500d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
7514b9ad928SBarry Smith /*@
7520d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
7534b9ad928SBarry Smith    complicated cycling.
7544b9ad928SBarry Smith 
755ad4df100SBarry Smith    Logically Collective on PC
7564b9ad928SBarry Smith 
7574b9ad928SBarry Smith    Input Parameters:
758c2be2410SBarry Smith +  pc - the multigrid context
7590d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
7604b9ad928SBarry Smith 
7614b9ad928SBarry Smith    Options Database Key:
7620d353602SBarry Smith $  -pc_mg_cycle_type v or w
7634b9ad928SBarry Smith 
7644b9ad928SBarry Smith    Level: advanced
7654b9ad928SBarry Smith 
7664b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7674b9ad928SBarry Smith 
7680d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
7694b9ad928SBarry Smith @*/
7707087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
7714b9ad928SBarry Smith {
772f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
773f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
77479416396SBarry Smith   PetscInt     i,levels;
7754b9ad928SBarry Smith 
7764b9ad928SBarry Smith   PetscFunctionBegin;
7770700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
778e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
779c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
780f3fbd535SBarry Smith   levels = mglevels[0]->levels;
7814b9ad928SBarry Smith 
7824b9ad928SBarry Smith   for (i=0; i<levels; i++) {
783f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
7844b9ad928SBarry Smith   }
7854b9ad928SBarry Smith   PetscFunctionReturn(0);
7864b9ad928SBarry Smith }
7874b9ad928SBarry Smith 
7884b9ad928SBarry Smith #undef __FUNCT__
7898cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
7908cc2d5dfSBarry Smith /*@
7918cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
7928cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
7938cc2d5dfSBarry Smith 
794ad4df100SBarry Smith    Logically Collective on PC
7958cc2d5dfSBarry Smith 
7968cc2d5dfSBarry Smith    Input Parameters:
7978cc2d5dfSBarry Smith +  pc - the multigrid context
7988cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
7998cc2d5dfSBarry Smith 
8008cc2d5dfSBarry Smith    Options Database Key:
8018cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
8028cc2d5dfSBarry Smith 
8038cc2d5dfSBarry Smith    Level: advanced
8048cc2d5dfSBarry Smith 
8058cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
8068cc2d5dfSBarry Smith 
8078cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8088cc2d5dfSBarry Smith 
8098cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
8108cc2d5dfSBarry Smith @*/
8117087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
8128cc2d5dfSBarry Smith {
813f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
814f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
8158cc2d5dfSBarry Smith   PetscInt     i,levels;
8168cc2d5dfSBarry Smith 
8178cc2d5dfSBarry Smith   PetscFunctionBegin;
8180700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
819e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
820c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
821f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8228cc2d5dfSBarry Smith 
8238cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
82431567311SBarry Smith     mg->cyclesperpcapply  = n;
8258cc2d5dfSBarry Smith   }
8268cc2d5dfSBarry Smith   PetscFunctionReturn(0);
8278cc2d5dfSBarry Smith }
8288cc2d5dfSBarry Smith 
8298cc2d5dfSBarry Smith #undef __FUNCT__
8309dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
831c2be2410SBarry Smith /*@
83297177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
833c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
834c2be2410SBarry Smith 
835ad4df100SBarry Smith    Logically Collective on PC
836c2be2410SBarry Smith 
837c2be2410SBarry Smith    Input Parameters:
8383fc8bf9cSBarry Smith .  pc - the multigrid context
839c2be2410SBarry Smith 
840c2be2410SBarry Smith    Options Database Key:
841c2be2410SBarry Smith $  -pc_mg_galerkin
842c2be2410SBarry Smith 
843c2be2410SBarry Smith    Level: intermediate
844c2be2410SBarry Smith 
845c2be2410SBarry Smith .keywords: MG, set, Galerkin
846c2be2410SBarry Smith 
8473fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
8483fc8bf9cSBarry Smith 
849c2be2410SBarry Smith @*/
8507087cfbeSBarry Smith PetscErrorCode  PCMGSetGalerkin(PC pc)
851c2be2410SBarry Smith {
852f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
853c2be2410SBarry Smith 
854c2be2410SBarry Smith   PetscFunctionBegin;
8550700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
856218a07d4SBarry Smith   mg->galerkin = PETSC_TRUE;
857c2be2410SBarry Smith   PetscFunctionReturn(0);
858c2be2410SBarry Smith }
859c2be2410SBarry Smith 
860c2be2410SBarry Smith #undef __FUNCT__
8613fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
8623fc8bf9cSBarry Smith /*@
8633fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
8643fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
8653fc8bf9cSBarry Smith 
8663fc8bf9cSBarry Smith    Not Collective
8673fc8bf9cSBarry Smith 
8683fc8bf9cSBarry Smith    Input Parameter:
8693fc8bf9cSBarry Smith .  pc - the multigrid context
8703fc8bf9cSBarry Smith 
8713fc8bf9cSBarry Smith    Output Parameter:
8723fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
8733fc8bf9cSBarry Smith 
8743fc8bf9cSBarry Smith    Options Database Key:
8753fc8bf9cSBarry Smith $  -pc_mg_galerkin
8763fc8bf9cSBarry Smith 
8773fc8bf9cSBarry Smith    Level: intermediate
8783fc8bf9cSBarry Smith 
8793fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
8803fc8bf9cSBarry Smith 
8813fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
8823fc8bf9cSBarry Smith 
8833fc8bf9cSBarry Smith @*/
8847087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
8853fc8bf9cSBarry Smith {
886f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
8873fc8bf9cSBarry Smith 
8883fc8bf9cSBarry Smith   PetscFunctionBegin;
8890700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
890218a07d4SBarry Smith   *galerkin = mg->galerkin;
8913fc8bf9cSBarry Smith   PetscFunctionReturn(0);
8923fc8bf9cSBarry Smith }
8933fc8bf9cSBarry Smith 
8943fc8bf9cSBarry Smith #undef __FUNCT__
8959dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
8964b9ad928SBarry Smith /*@
89797177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
89897177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
8994b9ad928SBarry Smith    pre-smoothing steps on different levels.
9004b9ad928SBarry Smith 
901ad4df100SBarry Smith    Logically Collective on PC
9024b9ad928SBarry Smith 
9034b9ad928SBarry Smith    Input Parameters:
9044b9ad928SBarry Smith +  mg - the multigrid context
9054b9ad928SBarry Smith -  n - the number of smoothing steps
9064b9ad928SBarry Smith 
9074b9ad928SBarry Smith    Options Database Key:
9084b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
9094b9ad928SBarry Smith 
9104b9ad928SBarry Smith    Level: advanced
9114b9ad928SBarry Smith 
9124b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
9134b9ad928SBarry Smith 
91497177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
9154b9ad928SBarry Smith @*/
9167087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
9174b9ad928SBarry Smith {
918f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
919f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9206849ba73SBarry Smith   PetscErrorCode ierr;
92179416396SBarry Smith   PetscInt       i,levels;
9224b9ad928SBarry Smith 
9234b9ad928SBarry Smith   PetscFunctionBegin;
9240700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
925e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
926c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
927f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9284b9ad928SBarry Smith 
929b05257ddSBarry Smith   for (i=1; i<levels; i++) {
9304b9ad928SBarry Smith     /* make sure smoother up and down are different */
93197177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
932f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
93331567311SBarry Smith     mg->default_smoothd = n;
9344b9ad928SBarry Smith   }
9354b9ad928SBarry Smith   PetscFunctionReturn(0);
9364b9ad928SBarry Smith }
9374b9ad928SBarry Smith 
9384b9ad928SBarry Smith #undef __FUNCT__
9399dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
9404b9ad928SBarry Smith /*@
94197177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
94297177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
9434b9ad928SBarry Smith    post-smoothing steps on different levels.
9444b9ad928SBarry Smith 
945ad4df100SBarry Smith    Logically Collective on PC
9464b9ad928SBarry Smith 
9474b9ad928SBarry Smith    Input Parameters:
9484b9ad928SBarry Smith +  mg - the multigrid context
9494b9ad928SBarry Smith -  n - the number of smoothing steps
9504b9ad928SBarry Smith 
9514b9ad928SBarry Smith    Options Database Key:
9524b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
9534b9ad928SBarry Smith 
9544b9ad928SBarry Smith    Level: advanced
9554b9ad928SBarry Smith 
9564b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
957a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
9584b9ad928SBarry Smith 
9594b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
9604b9ad928SBarry Smith 
96197177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
9624b9ad928SBarry Smith @*/
9637087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
9644b9ad928SBarry Smith {
965f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
966f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9676849ba73SBarry Smith   PetscErrorCode ierr;
96879416396SBarry Smith   PetscInt       i,levels;
9694b9ad928SBarry Smith 
9704b9ad928SBarry Smith   PetscFunctionBegin;
9710700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
972e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
973c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
974f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9754b9ad928SBarry Smith 
9764b9ad928SBarry Smith   for (i=1; i<levels; i++) {
9774b9ad928SBarry Smith     /* make sure smoother up and down are different */
97897177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
979f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
98031567311SBarry Smith     mg->default_smoothu = n;
9814b9ad928SBarry Smith   }
9824b9ad928SBarry Smith   PetscFunctionReturn(0);
9834b9ad928SBarry Smith }
9844b9ad928SBarry Smith 
9854b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
9864b9ad928SBarry Smith 
9873b09bd56SBarry Smith /*MC
988ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
9893b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
9903b09bd56SBarry Smith 
9913b09bd56SBarry Smith    Options Database Keys:
9923b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
9930d353602SBarry Smith .  -pc_mg_cycles v or w
99479416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
9953b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
9963b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
9973b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
9983b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
99968eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
10003b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1001e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
10023b09bd56SBarry Smith 
100324c3aa18SBarry 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
10043b09bd56SBarry Smith 
10053b09bd56SBarry Smith    Level: intermediate
10063b09bd56SBarry Smith 
10078f87f92bSBarry Smith    Concepts: multigrid/multilevel
10083b09bd56SBarry Smith 
100924c3aa18SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC
10100d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
101197177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
101297177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
10130d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
10143b09bd56SBarry Smith M*/
10153b09bd56SBarry Smith 
10164b9ad928SBarry Smith EXTERN_C_BEGIN
10174b9ad928SBarry Smith #undef __FUNCT__
10184b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
10197087cfbeSBarry Smith PetscErrorCode  PCCreate_MG(PC pc)
10204b9ad928SBarry Smith {
1021f3fbd535SBarry Smith   PC_MG          *mg;
1022f3fbd535SBarry Smith   PetscErrorCode ierr;
1023f3fbd535SBarry Smith 
10244b9ad928SBarry Smith   PetscFunctionBegin;
1025f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1026f3fbd535SBarry Smith   pc->data    = (void*)mg;
1027f3fbd535SBarry Smith   mg->nlevels = -1;
1028f3fbd535SBarry Smith 
10294b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
10304b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1031a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
10324b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
10334b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
10344b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
10354b9ad928SBarry Smith   PetscFunctionReturn(0);
10364b9ad928SBarry Smith }
10374b9ad928SBarry Smith EXTERN_C_END
1038