xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 2d2b81a6488da149505f49f1cb76c0ad7d19ae12)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
54b9ad928SBarry Smith */
67c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
74b9ad928SBarry Smith 
84b9ad928SBarry Smith 
94b9ad928SBarry Smith #undef __FUNCT__
109dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private"
1131567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
124b9ad928SBarry Smith {
1331567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1431567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
156849ba73SBarry Smith   PetscErrorCode ierr;
1631567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
174b9ad928SBarry Smith 
184b9ad928SBarry Smith   PetscFunctionBegin;
194b9ad928SBarry Smith 
2032cf1786SBarry Smith   if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2131567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
2232cf1786SBarry Smith   if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2331567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2432cf1786SBarry Smith     if (mg->eventresidual) {ierr = PetscLogEventBegin(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2531567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2632cf1786SBarry Smith     if (mg->eventresidual) {ierr = PetscLogEventEnd(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
274b9ad928SBarry Smith 
284b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
2931567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
304b9ad928SBarry Smith       PetscReal rnorm;
3131567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
324b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3370441072SBarry Smith         if (rnorm < mg->abstol) {
344d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
351e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
364b9ad928SBarry Smith         } else {
374d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
381e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr);
394b9ad928SBarry Smith         }
404b9ad928SBarry Smith         PetscFunctionReturn(0);
414b9ad928SBarry Smith       }
424b9ad928SBarry Smith     }
434b9ad928SBarry Smith 
4431567311SBarry Smith     mgc = *(mglevelsin - 1);
4532cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4631567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
4732cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
48efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
494b9ad928SBarry Smith     while (cycles--) {
5031567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
514b9ad928SBarry Smith     }
5232cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5331567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5432cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5532cf1786SBarry Smith     if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5631567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
5732cf1786SBarry Smith     if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
584b9ad928SBarry Smith   }
594b9ad928SBarry Smith   PetscFunctionReturn(0);
604b9ad928SBarry Smith }
614b9ad928SBarry Smith 
624b9ad928SBarry Smith #undef __FUNCT__
634b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
647319c654SBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
654b9ad928SBarry Smith {
66f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
67f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
68dfbe8321SBarry Smith   PetscErrorCode ierr;
69f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
704b9ad928SBarry Smith 
714b9ad928SBarry Smith   PetscFunctionBegin;
72f3fbd535SBarry Smith   mglevels[levels-1]->b    = b;
73f3fbd535SBarry Smith   mglevels[levels-1]->x    = x;
744b9ad928SBarry Smith 
7531567311SBarry Smith   mg->rtol = rtol;
7631567311SBarry Smith   mg->abstol = abstol;
7731567311SBarry Smith   mg->dtol = dtol;
784b9ad928SBarry Smith   if (rtol) {
794b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
804b9ad928SBarry Smith     PetscReal rnorm;
817319c654SBarry Smith     if (zeroguess) {
827319c654SBarry Smith       ierr               = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
837319c654SBarry Smith     } else {
84f3fbd535SBarry Smith       ierr               = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
854b9ad928SBarry Smith       ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
867319c654SBarry Smith     }
8731567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
8870441072SBarry Smith   } else if (abstol) {
8931567311SBarry Smith     mg->ttol = abstol;
904b9ad928SBarry Smith   } else {
9131567311SBarry Smith     mg->ttol = 0.0;
924b9ad928SBarry Smith   }
934b9ad928SBarry Smith 
944d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
954d0a8057SBarry Smith      stop prematurely do to small residual */
964d0a8057SBarry Smith   for (i=1; i<levels; i++) {
97f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
98f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
99f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1004b9ad928SBarry Smith     }
1014d0a8057SBarry Smith   }
1024d0a8057SBarry Smith 
1034d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1044d0a8057SBarry Smith   for (i=0; i<its; i++) {
105f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1064d0a8057SBarry Smith     if (*reason) break;
1074d0a8057SBarry Smith   }
1084d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1094d0a8057SBarry Smith   *outits = i;
1104b9ad928SBarry Smith   PetscFunctionReturn(0);
1114b9ad928SBarry Smith }
1124b9ad928SBarry Smith 
1134b9ad928SBarry Smith #undef __FUNCT__
1149dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
1154b9ad928SBarry Smith /*@C
11697177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1174b9ad928SBarry Smith    Must be called before any other MG routine.
1184b9ad928SBarry Smith 
1194b9ad928SBarry Smith    Collective on PC
1204b9ad928SBarry Smith 
1214b9ad928SBarry Smith    Input Parameters:
1224b9ad928SBarry Smith +  pc - the preconditioner context
1234b9ad928SBarry Smith .  levels - the number of levels
1244b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1254b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
1264b9ad928SBarry Smith 
1274b9ad928SBarry Smith    Level: intermediate
1284b9ad928SBarry Smith 
1294b9ad928SBarry Smith    Notes:
1304b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1314b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1324b9ad928SBarry Smith 
1334b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1344b9ad928SBarry Smith 
13597177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1364b9ad928SBarry Smith @*/
13797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1384b9ad928SBarry Smith {
139dfbe8321SBarry Smith   PetscErrorCode ierr;
140f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
141f3fbd535SBarry Smith   MPI_Comm       comm = ((PetscObject)pc)->comm;
142f3fbd535SBarry Smith   PC_MG_Levels   **mglevels;
143f3fbd535SBarry Smith   PetscInt       i;
144f3fbd535SBarry Smith   PetscMPIInt    size;
145f3fbd535SBarry Smith   const char     *prefix;
146f3fbd535SBarry Smith   PC             ipc;
1474b9ad928SBarry Smith 
1484b9ad928SBarry Smith   PetscFunctionBegin;
1490700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
150f3fbd535SBarry Smith   if (mg->nlevels > -1) {
151f3fbd535SBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n  make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
1524b9ad928SBarry Smith   }
153f3fbd535SBarry Smith   if (mg->levels) SETERRQ(PETSC_ERR_PLIB,"Internal error in PETSc, this array should not yet exist");
154f3fbd535SBarry Smith 
155f3fbd535SBarry Smith   mg->nlevels      = levels;
156218a07d4SBarry Smith   mg->galerkin     = PETSC_FALSE;
157218a07d4SBarry Smith   mg->galerkinused = PETSC_FALSE;
158f3fbd535SBarry Smith 
159f3fbd535SBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr);
160f3fbd535SBarry Smith   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
161f3fbd535SBarry Smith 
162f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
163f3fbd535SBarry Smith 
164f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
165f3fbd535SBarry Smith     ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr);
166f3fbd535SBarry Smith     mglevels[i]->level           = i;
167f3fbd535SBarry Smith     mglevels[i]->levels          = levels;
168f3fbd535SBarry Smith     mglevels[i]->cycles          = PC_MG_CYCLE_V;
16931567311SBarry Smith     mg->default_smoothu = 1;
17031567311SBarry Smith     mg->default_smoothd = 1;
171f3fbd535SBarry Smith 
172f3fbd535SBarry Smith     if (comms) comm = comms[i];
173f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
174f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
17531567311SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
176f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
177f3fbd535SBarry Smith 
178f3fbd535SBarry Smith     /* do special stuff for coarse grid */
179f3fbd535SBarry Smith     if (!i && levels > 1) {
180f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
181f3fbd535SBarry Smith 
182f3fbd535SBarry Smith       /* coarse solve is (redundant) LU by default */
183f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
184f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
185f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
186f3fbd535SBarry Smith       if (size > 1) {
187f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
188f3fbd535SBarry Smith       } else {
189f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
190f3fbd535SBarry Smith       }
191f3fbd535SBarry Smith 
192f3fbd535SBarry Smith     } else {
193f3fbd535SBarry Smith       char tprefix[128];
194f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
195f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
196f3fbd535SBarry Smith     }
197f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr);
198f3fbd535SBarry Smith     mglevels[i]->smoothu    = mglevels[i]->smoothd;
19931567311SBarry Smith     mg->rtol                = 0.0;
20031567311SBarry Smith     mg->abstol              = 0.0;
20131567311SBarry Smith     mg->dtol                = 0.0;
20231567311SBarry Smith     mg->ttol                = 0.0;
20331567311SBarry Smith     mg->eventsmoothsetup    = 0;
20431567311SBarry Smith     mg->eventsmoothsolve    = 0;
20531567311SBarry Smith     mg->eventresidual       = 0;
20631567311SBarry Smith     mg->eventinterprestrict = 0;
20731567311SBarry Smith     mg->cyclesperpcapply    = 1;
208f3fbd535SBarry Smith   }
20931567311SBarry Smith   mg->am          = PC_MG_MULTIPLICATIVE;
210f3fbd535SBarry Smith   mg->levels      = mglevels;
2114b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
2124b9ad928SBarry Smith   PetscFunctionReturn(0);
2134b9ad928SBarry Smith }
2144b9ad928SBarry Smith 
2154b9ad928SBarry Smith #undef __FUNCT__
216c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG_Private"
217c07bf074SBarry Smith PetscErrorCode PCDestroy_MG_Private(PC pc)
218f3fbd535SBarry Smith {
219f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
220f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
221f3fbd535SBarry Smith   PetscErrorCode ierr;
222f3fbd535SBarry Smith   PetscInt       i,n;
223f3fbd535SBarry Smith 
224f3fbd535SBarry Smith   PetscFunctionBegin;
225f3fbd535SBarry Smith   if (mglevels) {
226f3fbd535SBarry Smith     n = mglevels[0]->levels;
227f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
228f3fbd535SBarry Smith       if (mglevels[i+1]->r) {ierr = VecDestroy(mglevels[i+1]->r);CHKERRQ(ierr);}
229f3fbd535SBarry Smith       if (mglevels[i]->b) {ierr = VecDestroy(mglevels[i]->b);CHKERRQ(ierr);}
230f3fbd535SBarry Smith       if (mglevels[i]->x) {ierr = VecDestroy(mglevels[i]->x);CHKERRQ(ierr);}
231f3fbd535SBarry Smith       if (mglevels[i+1]->restrct) {ierr = MatDestroy(mglevels[i+1]->restrct);CHKERRQ(ierr);}
232f3fbd535SBarry Smith       if (mglevels[i+1]->interpolate) {ierr = MatDestroy(mglevels[i+1]->interpolate);CHKERRQ(ierr);}
233f3fbd535SBarry Smith     }
234f3fbd535SBarry Smith 
235f3fbd535SBarry Smith     for (i=0; i<n; i++) {
236f3fbd535SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
237f3fbd535SBarry Smith 	ierr = KSPDestroy(mglevels[i]->smoothd);CHKERRQ(ierr);
238f3fbd535SBarry Smith       }
239f3fbd535SBarry Smith       ierr = KSPDestroy(mglevels[i]->smoothu);CHKERRQ(ierr);
240f3fbd535SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
241f3fbd535SBarry Smith     }
242f3fbd535SBarry Smith     ierr = PetscFree(mglevels);CHKERRQ(ierr);
243f3fbd535SBarry Smith   }
244c07bf074SBarry Smith   mg->nlevels = -1;
245c07bf074SBarry Smith   mg->levels  = PETSC_NULL;
246c07bf074SBarry Smith   PetscFunctionReturn(0);
247c07bf074SBarry Smith }
248c07bf074SBarry Smith 
249c07bf074SBarry Smith #undef __FUNCT__
250c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
251c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
252c07bf074SBarry Smith {
253c07bf074SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
254c07bf074SBarry Smith   PetscErrorCode ierr;
255c07bf074SBarry Smith 
256c07bf074SBarry Smith   PetscFunctionBegin;
257c07bf074SBarry Smith   ierr = PCDestroy_MG_Private(pc);CHKERRQ(ierr);
258f3fbd535SBarry Smith   ierr = PetscFree(mg);CHKERRQ(ierr);
259f3fbd535SBarry Smith   PetscFunctionReturn(0);
260f3fbd535SBarry Smith }
261f3fbd535SBarry Smith 
262f3fbd535SBarry Smith 
263f3fbd535SBarry Smith 
26431567311SBarry Smith EXTERN PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
265f3fbd535SBarry Smith EXTERN PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
26631567311SBarry Smith EXTERN PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
267f3fbd535SBarry Smith 
268f3fbd535SBarry Smith /*
269f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
270f3fbd535SBarry Smith              or full cycle of multigrid.
271f3fbd535SBarry Smith 
272f3fbd535SBarry Smith   Note:
273f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
274f3fbd535SBarry Smith */
275f3fbd535SBarry Smith #undef __FUNCT__
276f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
277f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
278f3fbd535SBarry Smith {
279f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
280f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
281f3fbd535SBarry Smith   PetscErrorCode ierr;
282f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
283f3fbd535SBarry Smith 
284f3fbd535SBarry Smith   PetscFunctionBegin;
285f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
286f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
28731567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
288f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
28931567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
290f3fbd535SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
291f3fbd535SBarry Smith     }
292f3fbd535SBarry Smith   }
29331567311SBarry Smith   else if (mg->am == PC_MG_ADDITIVE) {
29431567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
295f3fbd535SBarry Smith   }
29631567311SBarry Smith   else if (mg->am == PC_MG_KASKADE) {
29731567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
298f3fbd535SBarry Smith   }
299f3fbd535SBarry Smith   else {
300f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
301f3fbd535SBarry Smith   }
302f3fbd535SBarry Smith   PetscFunctionReturn(0);
303f3fbd535SBarry Smith }
304f3fbd535SBarry Smith 
305f3fbd535SBarry Smith 
306f3fbd535SBarry Smith #undef __FUNCT__
307f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
308f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc)
309f3fbd535SBarry Smith {
310f3fbd535SBarry Smith   PetscErrorCode ierr;
311f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
312f3fbd535SBarry Smith   PetscTruth     flg;
313f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
314f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
315f3fbd535SBarry Smith   PCMGType       mgtype;
316f3fbd535SBarry Smith   PCMGCycleType  mgctype;
317f3fbd535SBarry Smith 
318f3fbd535SBarry Smith   PetscFunctionBegin;
319f3fbd535SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
32018aabeadSBarry Smith     if (!mglevels) {
321f3fbd535SBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
322f3fbd535SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
323f3fbd535SBarry Smith       mglevels = mg->levels;
324f3fbd535SBarry Smith     }
325f3fbd535SBarry Smith     mgctype = (PCMGCycleType) mglevels[0]->cycles;
326f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
327f3fbd535SBarry Smith     if (flg) {
328f3fbd535SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
329f3fbd535SBarry Smith     };
330f3fbd535SBarry Smith     flg  = PETSC_FALSE;
331f3fbd535SBarry Smith     ierr = PetscOptionsTruth("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
332f3fbd535SBarry Smith     if (flg) {
333f3fbd535SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
334f3fbd535SBarry Smith     }
335f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
336f3fbd535SBarry Smith     if (flg) {
337f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
338f3fbd535SBarry Smith     }
339f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
340f3fbd535SBarry Smith     if (flg) {
341f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
342f3fbd535SBarry Smith     }
34331567311SBarry Smith     mgtype = mg->am;
344f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
345f3fbd535SBarry Smith     if (flg) {
346f3fbd535SBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
347f3fbd535SBarry Smith     }
34831567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
34931567311SBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
350f3fbd535SBarry Smith       if (flg) {
351f3fbd535SBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
352f3fbd535SBarry Smith       }
353f3fbd535SBarry Smith     }
354f3fbd535SBarry Smith     flg  = PETSC_FALSE;
355f3fbd535SBarry Smith     ierr = PetscOptionsTruth("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
356f3fbd535SBarry Smith     if (flg) {
357f3fbd535SBarry Smith       PetscInt i;
358f3fbd535SBarry Smith       char     eventname[128];
359f3fbd535SBarry Smith       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
360f3fbd535SBarry Smith       levels = mglevels[0]->levels;
361f3fbd535SBarry Smith       for (i=0; i<levels; i++) {
362f3fbd535SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
3630700a824SBarry Smith         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mg->eventsmoothsetup);CHKERRQ(ierr);
364f3fbd535SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
3650700a824SBarry Smith         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mg->eventsmoothsolve);CHKERRQ(ierr);
366f3fbd535SBarry Smith         if (i) {
367f3fbd535SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
3680700a824SBarry Smith           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mg->eventresidual);CHKERRQ(ierr);
369f3fbd535SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
3700700a824SBarry Smith           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mg->eventinterprestrict);CHKERRQ(ierr);
371f3fbd535SBarry Smith         }
372f3fbd535SBarry Smith       }
373f3fbd535SBarry Smith     }
374f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
375f3fbd535SBarry Smith   PetscFunctionReturn(0);
376f3fbd535SBarry Smith }
377f3fbd535SBarry Smith 
378f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
379f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
380f3fbd535SBarry Smith 
381f3fbd535SBarry Smith #undef __FUNCT__
382f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
383f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
384f3fbd535SBarry Smith {
385f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
386f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
387f3fbd535SBarry Smith   PetscErrorCode ierr;
388f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
389f3fbd535SBarry Smith   PetscTruth     iascii;
390f3fbd535SBarry Smith 
391f3fbd535SBarry Smith   PetscFunctionBegin;
392f3fbd535SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
393f3fbd535SBarry Smith   if (iascii) {
39431567311SBarry 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);
39531567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
39631567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
397f3fbd535SBarry Smith     }
398218a07d4SBarry Smith     if (mg->galerkin) {
399f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
400f3fbd535SBarry Smith     }
401f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
402f3fbd535SBarry Smith       if (!i) {
40331567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D presmooths=%D postsmooths=%D -----\n",i,mg->default_smoothd,mg->default_smoothu);CHKERRQ(ierr);
404f3fbd535SBarry Smith       } else {
40531567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
406f3fbd535SBarry Smith       }
407f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
408f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
409f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
410f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
411f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
412f3fbd535SBarry Smith       } else if (i){
41331567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr);
414f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
415f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
416f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
417f3fbd535SBarry Smith       }
418f3fbd535SBarry Smith     }
419f3fbd535SBarry Smith   } else {
420f3fbd535SBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
421f3fbd535SBarry Smith   }
422f3fbd535SBarry Smith   PetscFunctionReturn(0);
423f3fbd535SBarry Smith }
424f3fbd535SBarry Smith 
425f3fbd535SBarry Smith /*
426f3fbd535SBarry Smith     Calls setup for the KSP on each level
427f3fbd535SBarry Smith */
428f3fbd535SBarry Smith #undef __FUNCT__
429f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
430f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
431f3fbd535SBarry Smith {
432f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
433f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
434f3fbd535SBarry Smith   PetscErrorCode          ierr;
435f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
436f3fbd535SBarry Smith   PC                      cpc,mpc;
437f3fbd535SBarry Smith   PetscTruth              preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset;
438f3fbd535SBarry Smith   PetscViewerASCIIMonitor ascii;
439f3fbd535SBarry Smith   PetscViewer             viewer = PETSC_NULL;
440f3fbd535SBarry Smith   MPI_Comm                comm;
441f3fbd535SBarry Smith   Mat                     dA,dB;
442f3fbd535SBarry Smith   MatStructure            uflag;
443f3fbd535SBarry Smith   Vec                     tvec;
444218a07d4SBarry Smith   DM                      *dms;
445f3fbd535SBarry Smith 
446f3fbd535SBarry Smith   PetscFunctionBegin;
447f3fbd535SBarry Smith 
448f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
449f3fbd535SBarry Smith   /* so use those from global PC */
450f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
451f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
452f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
453f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr);
4541338a6b9SJed Brown   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) {
455f3fbd535SBarry 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);
456f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
457f3fbd535SBarry Smith   }
458f3fbd535SBarry Smith 
459*2d2b81a6SBarry Smith   if (pc->dm && !pc->setupcalled) {
460*2d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
461*2d2b81a6SBarry Smith     Mat A,p;
462218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
463218a07d4SBarry Smith     dms[n-1] = pc->dm;
464218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
465218a07d4SBarry Smith       ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr);
466*2d2b81a6SBarry Smith       ierr = DMGetInterpolation(dms[i],dms[i+1],&p,PETSC_NULL);CHKERRQ(ierr);
467*2d2b81a6SBarry Smith       ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
468*2d2b81a6SBarry Smith       ierr = MatDestroy(p);CHKERRQ(ierr);
469218a07d4SBarry Smith     }
470*2d2b81a6SBarry Smith 
471*2d2b81a6SBarry Smith     if (!mg->galerkin) {
472*2d2b81a6SBarry Smith       /* build the empty matrices for the coarser meshes */
473*2d2b81a6SBarry Smith       for (i=n-2; i>-1; i--) {
474*2d2b81a6SBarry Smith 	ierr = DMGetMatrix(dms[i],MATAIJ,&A);CHKERRQ(ierr);
475*2d2b81a6SBarry Smith 	ierr = KSPSetOperators(mglevels[i]->smoothd,A,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
476*2d2b81a6SBarry Smith 	ierr = MatDestroy(A);CHKERRQ(ierr);
477*2d2b81a6SBarry Smith       }
478*2d2b81a6SBarry Smith     }
479*2d2b81a6SBarry Smith 
480218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
481218a07d4SBarry Smith       ierr = DMDestroy(dms[i]);CHKERRQ(ierr);
482218a07d4SBarry Smith     }
483*2d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
484218a07d4SBarry Smith   }
485218a07d4SBarry Smith 
486218a07d4SBarry Smith   if (mg->galerkin) {
487f3fbd535SBarry Smith     Mat B;
488218a07d4SBarry Smith     mg->galerkinused = PETSC_TRUE;
489f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
490f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
491f3fbd535SBarry Smith     if (!pc->setupcalled) {
492f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
493f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
494f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
495f3fbd535SBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
496f3fbd535SBarry Smith         dB   = B;
497f3fbd535SBarry Smith       }
498f3fbd535SBarry Smith       ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);
499f3fbd535SBarry Smith     } else {
500f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
501f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
502f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
503f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
504f3fbd535SBarry Smith         dB   = B;
505f3fbd535SBarry Smith       }
506f3fbd535SBarry Smith     }
507f3fbd535SBarry Smith   }
508f3fbd535SBarry Smith 
509f3fbd535SBarry Smith   if (!pc->setupcalled) {
510f3fbd535SBarry Smith     ierr = PetscOptionsGetTruth(0,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr);
511f3fbd535SBarry Smith 
512f3fbd535SBarry Smith     for (i=0; i<n; i++) {
513f3fbd535SBarry Smith       if (monitor) {
514f3fbd535SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr);
515f3fbd535SBarry Smith         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
516f3fbd535SBarry Smith         ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
517f3fbd535SBarry Smith       }
518f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
519f3fbd535SBarry Smith     }
520f3fbd535SBarry Smith     for (i=1; i<n; i++) {
521f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
522f3fbd535SBarry Smith         if (monitor) {
523f3fbd535SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr);
524f3fbd535SBarry Smith           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
525f3fbd535SBarry Smith           ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
526f3fbd535SBarry Smith         }
527f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
528f3fbd535SBarry Smith       }
529f3fbd535SBarry Smith     }
530f3fbd535SBarry Smith     for (i=1; i<n; i++) {
531f3fbd535SBarry Smith       if (!mglevels[i]->residual) {
532f3fbd535SBarry Smith         Mat mat;
533f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
534f3fbd535SBarry Smith         ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
535f3fbd535SBarry Smith       }
536f3fbd535SBarry Smith       if (mglevels[i]->restrct && !mglevels[i]->interpolate) {
537f3fbd535SBarry Smith         ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr);
538f3fbd535SBarry Smith       }
539f3fbd535SBarry Smith       if (!mglevels[i]->restrct && mglevels[i]->interpolate) {
540f3fbd535SBarry Smith         ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr);
541f3fbd535SBarry Smith       }
542f3fbd535SBarry Smith #if defined(PETSC_USE_DEBUG)
543f3fbd535SBarry Smith       if (!mglevels[i]->restrct || !mglevels[i]->interpolate) {
544f3fbd535SBarry Smith         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
545f3fbd535SBarry Smith       }
546f3fbd535SBarry Smith #endif
547f3fbd535SBarry Smith     }
548f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
549f3fbd535SBarry Smith       if (!mglevels[i]->b) {
550f3fbd535SBarry Smith         Vec *vec;
551f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
552f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
553f3fbd535SBarry Smith         ierr = VecDestroy(*vec);CHKERRQ(ierr);
554f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
555f3fbd535SBarry Smith       }
556f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
557f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
558f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
559f3fbd535SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
560f3fbd535SBarry Smith       }
561f3fbd535SBarry Smith       if (!mglevels[i]->x) {
562f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
563f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
564f3fbd535SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
565f3fbd535SBarry Smith       }
566f3fbd535SBarry Smith     }
567f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
568f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
569f3fbd535SBarry Smith       Vec *vec;
570f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
571f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
572f3fbd535SBarry Smith       ierr = VecDestroy(*vec);CHKERRQ(ierr);
573f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
574f3fbd535SBarry Smith     }
575f3fbd535SBarry Smith   }
576f3fbd535SBarry Smith 
577f3fbd535SBarry Smith 
578f3fbd535SBarry Smith   for (i=1; i<n; i++) {
579f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
580f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
581f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
582f3fbd535SBarry Smith     }
58331567311SBarry Smith     if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
584f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
58531567311SBarry Smith     if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
586f3fbd535SBarry Smith   }
587f3fbd535SBarry Smith   for (i=1; i<n; i++) {
588f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
589f3fbd535SBarry Smith       Mat          downmat,downpmat;
590f3fbd535SBarry Smith       MatStructure matflag;
591f3fbd535SBarry Smith       PetscTruth   opsset;
592f3fbd535SBarry Smith 
593f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
594f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
595f3fbd535SBarry Smith       if (!opsset) {
596f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
597f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
598f3fbd535SBarry Smith       }
599f3fbd535SBarry Smith 
600f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
60131567311SBarry Smith       if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
602f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
60331567311SBarry Smith       if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
604f3fbd535SBarry Smith     }
605f3fbd535SBarry Smith   }
606f3fbd535SBarry Smith 
607f3fbd535SBarry Smith   /*
608f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
609f3fbd535SBarry Smith   */
610f3fbd535SBarry Smith   ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
611f3fbd535SBarry Smith   if (preonly) {
612f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
613f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
614f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
615f3fbd535SBarry Smith     if (!lu && !redundant && !cholesky) {
616f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
617f3fbd535SBarry Smith     }
618f3fbd535SBarry Smith   }
619f3fbd535SBarry Smith 
620f3fbd535SBarry Smith   if (!pc->setupcalled) {
621f3fbd535SBarry Smith     if (monitor) {
622f3fbd535SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr);
623f3fbd535SBarry Smith       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
624f3fbd535SBarry Smith       ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
625f3fbd535SBarry Smith     }
626f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
627f3fbd535SBarry Smith   }
628f3fbd535SBarry Smith 
62931567311SBarry Smith   if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
630f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
63131567311SBarry Smith   if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
632f3fbd535SBarry Smith 
633f3fbd535SBarry Smith   /*
634f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
635f3fbd535SBarry Smith    Jacobian/stiffness on each level. This allows Matlab users to
636f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
637f3fbd535SBarry Smith 
638f3fbd535SBarry Smith    Only support one or the other at the same time.
639f3fbd535SBarry Smith   */
640f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
641f3fbd535SBarry Smith   ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
642f3fbd535SBarry Smith   if (dump) {
643f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
644f3fbd535SBarry Smith   }
645f3fbd535SBarry Smith   dump = PETSC_FALSE;
646f3fbd535SBarry Smith #endif
647f3fbd535SBarry Smith   ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
648f3fbd535SBarry Smith   if (dump) {
649f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
650f3fbd535SBarry Smith   }
651f3fbd535SBarry Smith 
652f3fbd535SBarry Smith   if (viewer) {
653f3fbd535SBarry Smith     for (i=1; i<n; i++) {
654f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
655f3fbd535SBarry Smith     }
656f3fbd535SBarry Smith     for (i=0; i<n; i++) {
657f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
658f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
659f3fbd535SBarry Smith     }
660f3fbd535SBarry Smith   }
661f3fbd535SBarry Smith   PetscFunctionReturn(0);
662f3fbd535SBarry Smith }
663f3fbd535SBarry Smith 
664f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
665f3fbd535SBarry Smith 
666f3fbd535SBarry Smith #undef __FUNCT__
6679dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
6684b9ad928SBarry Smith /*@
66997177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
6704b9ad928SBarry Smith 
6714b9ad928SBarry Smith    Not Collective
6724b9ad928SBarry Smith 
6734b9ad928SBarry Smith    Input Parameter:
6744b9ad928SBarry Smith .  pc - the preconditioner context
6754b9ad928SBarry Smith 
6764b9ad928SBarry Smith    Output parameter:
6774b9ad928SBarry Smith .  levels - the number of levels
6784b9ad928SBarry Smith 
6794b9ad928SBarry Smith    Level: advanced
6804b9ad928SBarry Smith 
6814b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
6824b9ad928SBarry Smith 
68397177400SBarry Smith .seealso: PCMGSetLevels()
6844b9ad928SBarry Smith @*/
68597177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
6864b9ad928SBarry Smith {
687f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
6884b9ad928SBarry Smith 
6894b9ad928SBarry Smith   PetscFunctionBegin;
6900700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
6914482741eSBarry Smith   PetscValidIntPointer(levels,2);
692f3fbd535SBarry Smith   *levels = mg->nlevels;
6934b9ad928SBarry Smith   PetscFunctionReturn(0);
6944b9ad928SBarry Smith }
6954b9ad928SBarry Smith 
6964b9ad928SBarry Smith #undef __FUNCT__
6979dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
6984b9ad928SBarry Smith /*@
69997177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
7004b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
7014b9ad928SBarry Smith 
7024b9ad928SBarry Smith    Collective on PC
7034b9ad928SBarry Smith 
7044b9ad928SBarry Smith    Input Parameters:
7054b9ad928SBarry Smith +  pc - the preconditioner context
7069dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
7079dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
7084b9ad928SBarry Smith 
7094b9ad928SBarry Smith    Options Database Key:
7104b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
7114b9ad928SBarry Smith    additive, full, kaskade
7124b9ad928SBarry Smith 
7134b9ad928SBarry Smith    Level: advanced
7144b9ad928SBarry Smith 
7154b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
7164b9ad928SBarry Smith 
71797177400SBarry Smith .seealso: PCMGSetLevels()
7184b9ad928SBarry Smith @*/
7199dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
7204b9ad928SBarry Smith {
721f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
7224b9ad928SBarry Smith 
7234b9ad928SBarry Smith   PetscFunctionBegin;
7240700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
72531567311SBarry Smith   mg->am = form;
7269dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
7274b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
7284b9ad928SBarry Smith   PetscFunctionReturn(0);
7294b9ad928SBarry Smith }
7304b9ad928SBarry Smith 
7314b9ad928SBarry Smith #undef __FUNCT__
7320d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
7334b9ad928SBarry Smith /*@
7340d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
7354b9ad928SBarry Smith    complicated cycling.
7364b9ad928SBarry Smith 
7374b9ad928SBarry Smith    Collective on PC
7384b9ad928SBarry Smith 
7394b9ad928SBarry Smith    Input Parameters:
740c2be2410SBarry Smith +  pc - the multigrid context
7410d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
7424b9ad928SBarry Smith 
7434b9ad928SBarry Smith    Options Database Key:
7440d353602SBarry Smith $  -pc_mg_cycle_type v or w
7454b9ad928SBarry Smith 
7464b9ad928SBarry Smith    Level: advanced
7474b9ad928SBarry Smith 
7484b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7494b9ad928SBarry Smith 
7500d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
7514b9ad928SBarry Smith @*/
7520d353602SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n)
7534b9ad928SBarry Smith {
754f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
755f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
75679416396SBarry Smith   PetscInt     i,levels;
7574b9ad928SBarry Smith 
7584b9ad928SBarry Smith   PetscFunctionBegin;
7590700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
760f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
761f3fbd535SBarry Smith   levels = mglevels[0]->levels;
7624b9ad928SBarry Smith 
7634b9ad928SBarry Smith   for (i=0; i<levels; i++) {
764f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
7654b9ad928SBarry Smith   }
7664b9ad928SBarry Smith   PetscFunctionReturn(0);
7674b9ad928SBarry Smith }
7684b9ad928SBarry Smith 
7694b9ad928SBarry Smith #undef __FUNCT__
7708cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
7718cc2d5dfSBarry Smith /*@
7728cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
7738cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
7748cc2d5dfSBarry Smith 
7758cc2d5dfSBarry Smith    Collective on PC
7768cc2d5dfSBarry Smith 
7778cc2d5dfSBarry Smith    Input Parameters:
7788cc2d5dfSBarry Smith +  pc - the multigrid context
7798cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
7808cc2d5dfSBarry Smith 
7818cc2d5dfSBarry Smith    Options Database Key:
7828cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
7838cc2d5dfSBarry Smith 
7848cc2d5dfSBarry Smith    Level: advanced
7858cc2d5dfSBarry Smith 
7868cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
7878cc2d5dfSBarry Smith 
7888cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7898cc2d5dfSBarry Smith 
7908cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
7918cc2d5dfSBarry Smith @*/
7928cc2d5dfSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
7938cc2d5dfSBarry Smith {
794f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
795f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
7968cc2d5dfSBarry Smith   PetscInt     i,levels;
7978cc2d5dfSBarry Smith 
7988cc2d5dfSBarry Smith   PetscFunctionBegin;
7990700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
800f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
801f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8028cc2d5dfSBarry Smith 
8038cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
80431567311SBarry Smith     mg->cyclesperpcapply  = n;
8058cc2d5dfSBarry Smith   }
8068cc2d5dfSBarry Smith   PetscFunctionReturn(0);
8078cc2d5dfSBarry Smith }
8088cc2d5dfSBarry Smith 
8098cc2d5dfSBarry Smith #undef __FUNCT__
8109dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
811c2be2410SBarry Smith /*@
81297177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
813c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
814c2be2410SBarry Smith 
815c2be2410SBarry Smith    Collective on PC
816c2be2410SBarry Smith 
817c2be2410SBarry Smith    Input Parameters:
8183fc8bf9cSBarry Smith .  pc - the multigrid context
819c2be2410SBarry Smith 
820c2be2410SBarry Smith    Options Database Key:
821c2be2410SBarry Smith $  -pc_mg_galerkin
822c2be2410SBarry Smith 
823c2be2410SBarry Smith    Level: intermediate
824c2be2410SBarry Smith 
825c2be2410SBarry Smith .keywords: MG, set, Galerkin
826c2be2410SBarry Smith 
8273fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
8283fc8bf9cSBarry Smith 
829c2be2410SBarry Smith @*/
83097177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
831c2be2410SBarry Smith {
832f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
833c2be2410SBarry Smith 
834c2be2410SBarry Smith   PetscFunctionBegin;
8350700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
836218a07d4SBarry Smith   mg->galerkin = PETSC_TRUE;
837c2be2410SBarry Smith   PetscFunctionReturn(0);
838c2be2410SBarry Smith }
839c2be2410SBarry Smith 
840c2be2410SBarry Smith #undef __FUNCT__
8413fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
8423fc8bf9cSBarry Smith /*@
8433fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
8443fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
8453fc8bf9cSBarry Smith 
8463fc8bf9cSBarry Smith    Not Collective
8473fc8bf9cSBarry Smith 
8483fc8bf9cSBarry Smith    Input Parameter:
8493fc8bf9cSBarry Smith .  pc - the multigrid context
8503fc8bf9cSBarry Smith 
8513fc8bf9cSBarry Smith    Output Parameter:
8523fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
8533fc8bf9cSBarry Smith 
8543fc8bf9cSBarry Smith    Options Database Key:
8553fc8bf9cSBarry Smith $  -pc_mg_galerkin
8563fc8bf9cSBarry Smith 
8573fc8bf9cSBarry Smith    Level: intermediate
8583fc8bf9cSBarry Smith 
8593fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
8603fc8bf9cSBarry Smith 
8613fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
8623fc8bf9cSBarry Smith 
8633fc8bf9cSBarry Smith @*/
8643fc8bf9cSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
8653fc8bf9cSBarry Smith {
866f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
8673fc8bf9cSBarry Smith 
8683fc8bf9cSBarry Smith   PetscFunctionBegin;
8690700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
870218a07d4SBarry Smith   *galerkin = mg->galerkin;
8713fc8bf9cSBarry Smith   PetscFunctionReturn(0);
8723fc8bf9cSBarry Smith }
8733fc8bf9cSBarry Smith 
8743fc8bf9cSBarry Smith #undef __FUNCT__
8759dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
8764b9ad928SBarry Smith /*@
87797177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
87897177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
8794b9ad928SBarry Smith    pre-smoothing steps on different levels.
8804b9ad928SBarry Smith 
8814b9ad928SBarry Smith    Collective on PC
8824b9ad928SBarry Smith 
8834b9ad928SBarry Smith    Input Parameters:
8844b9ad928SBarry Smith +  mg - the multigrid context
8854b9ad928SBarry Smith -  n - the number of smoothing steps
8864b9ad928SBarry Smith 
8874b9ad928SBarry Smith    Options Database Key:
8884b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
8894b9ad928SBarry Smith 
8904b9ad928SBarry Smith    Level: advanced
8914b9ad928SBarry Smith 
8924b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
8934b9ad928SBarry Smith 
89497177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
8954b9ad928SBarry Smith @*/
89697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
8974b9ad928SBarry Smith {
898f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
899f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9006849ba73SBarry Smith   PetscErrorCode ierr;
90179416396SBarry Smith   PetscInt       i,levels;
9024b9ad928SBarry Smith 
9034b9ad928SBarry Smith   PetscFunctionBegin;
9040700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
905f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
906f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9074b9ad928SBarry Smith 
908b05257ddSBarry Smith   for (i=1; i<levels; i++) {
9094b9ad928SBarry Smith     /* make sure smoother up and down are different */
91097177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
911f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
91231567311SBarry Smith     mg->default_smoothd = n;
9134b9ad928SBarry Smith   }
9144b9ad928SBarry Smith   PetscFunctionReturn(0);
9154b9ad928SBarry Smith }
9164b9ad928SBarry Smith 
9174b9ad928SBarry Smith #undef __FUNCT__
9189dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
9194b9ad928SBarry Smith /*@
92097177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
92197177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
9224b9ad928SBarry Smith    post-smoothing steps on different levels.
9234b9ad928SBarry Smith 
9244b9ad928SBarry Smith    Collective on PC
9254b9ad928SBarry Smith 
9264b9ad928SBarry Smith    Input Parameters:
9274b9ad928SBarry Smith +  mg - the multigrid context
9284b9ad928SBarry Smith -  n - the number of smoothing steps
9294b9ad928SBarry Smith 
9304b9ad928SBarry Smith    Options Database Key:
9314b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
9324b9ad928SBarry Smith 
9334b9ad928SBarry Smith    Level: advanced
9344b9ad928SBarry Smith 
9354b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
936a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
9374b9ad928SBarry Smith 
9384b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
9394b9ad928SBarry Smith 
94097177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
9414b9ad928SBarry Smith @*/
94297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
9434b9ad928SBarry Smith {
944f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
945f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9466849ba73SBarry Smith   PetscErrorCode ierr;
94779416396SBarry Smith   PetscInt       i,levels;
9484b9ad928SBarry Smith 
9494b9ad928SBarry Smith   PetscFunctionBegin;
9500700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
951f3fbd535SBarry Smith   if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
952f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9534b9ad928SBarry Smith 
9544b9ad928SBarry Smith   for (i=1; i<levels; i++) {
9554b9ad928SBarry Smith     /* make sure smoother up and down are different */
95697177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
957f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
95831567311SBarry Smith     mg->default_smoothu = n;
9594b9ad928SBarry Smith   }
9604b9ad928SBarry Smith   PetscFunctionReturn(0);
9614b9ad928SBarry Smith }
9624b9ad928SBarry Smith 
9634b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
9644b9ad928SBarry Smith 
9653b09bd56SBarry Smith /*MC
966ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
9673b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
9683b09bd56SBarry Smith 
9693b09bd56SBarry Smith    Options Database Keys:
9703b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
9710d353602SBarry Smith .  -pc_mg_cycles v or w
97279416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
9733b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
9743b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
9753b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
9763b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
97768eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
9783b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
9793b09bd56SBarry Smith                         to the Socket viewer for reading from Matlab.
9803b09bd56SBarry Smith 
9813b09bd56SBarry Smith    Notes:
9823b09bd56SBarry Smith 
9833b09bd56SBarry Smith    Level: intermediate
9843b09bd56SBarry Smith 
9858f87f92bSBarry Smith    Concepts: multigrid/multilevel
9863b09bd56SBarry Smith 
9873b09bd56SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
9880d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
98997177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
99097177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
9910d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
9923b09bd56SBarry Smith M*/
9933b09bd56SBarry Smith 
9944b9ad928SBarry Smith EXTERN_C_BEGIN
9954b9ad928SBarry Smith #undef __FUNCT__
9964b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
997dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
9984b9ad928SBarry Smith {
999f3fbd535SBarry Smith   PC_MG          *mg;
1000f3fbd535SBarry Smith   PetscErrorCode ierr;
1001f3fbd535SBarry Smith 
10024b9ad928SBarry Smith   PetscFunctionBegin;
1003f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1004f3fbd535SBarry Smith   pc->data    = (void*)mg;
1005f3fbd535SBarry Smith   mg->nlevels = -1;
1006f3fbd535SBarry Smith 
10074b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
10084b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
10094b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
10104b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
10114b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
10124b9ad928SBarry Smith   PetscFunctionReturn(0);
10134b9ad928SBarry Smith }
10144b9ad928SBarry Smith EXTERN_C_END
1015