xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 72a02f061ae258f536e9a337d12e6c620e12e349)
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
94336babb1SJed Brown      stop prematurely due 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__
1133aeaf226SBarry Smith #define __FUNCT__ "PCReset_MG"
1143aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1153aeaf226SBarry Smith {
1163aeaf226SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1173aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1183aeaf226SBarry Smith   PetscErrorCode ierr;
1193aeaf226SBarry Smith   PetscInt       i,n;
1203aeaf226SBarry Smith 
1213aeaf226SBarry Smith   PetscFunctionBegin;
1223aeaf226SBarry Smith   if (mglevels) {
1233aeaf226SBarry Smith     n = mglevels[0]->levels;
1243aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1253aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1263aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1273aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
1283aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1293aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
13073250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1313aeaf226SBarry Smith     }
1323aeaf226SBarry Smith 
1333aeaf226SBarry Smith     for (i=0; i<n; i++) {
1343aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1353aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1363aeaf226SBarry Smith 	ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1373aeaf226SBarry Smith       }
1383aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1393aeaf226SBarry Smith     }
1403aeaf226SBarry Smith   }
1413aeaf226SBarry Smith   PetscFunctionReturn(0);
1423aeaf226SBarry Smith }
1433aeaf226SBarry Smith 
1443aeaf226SBarry Smith #undef __FUNCT__
1459dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
1464b9ad928SBarry Smith /*@C
14797177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1484b9ad928SBarry Smith    Must be called before any other MG routine.
1494b9ad928SBarry Smith 
150ad4df100SBarry Smith    Logically Collective on PC
1514b9ad928SBarry Smith 
1524b9ad928SBarry Smith    Input Parameters:
1534b9ad928SBarry Smith +  pc - the preconditioner context
1544b9ad928SBarry Smith .  levels - the number of levels
1554b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1564b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
1574b9ad928SBarry Smith 
1584b9ad928SBarry Smith    Level: intermediate
1594b9ad928SBarry Smith 
1604b9ad928SBarry Smith    Notes:
1614b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1624b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1634b9ad928SBarry Smith 
1644b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1654b9ad928SBarry Smith 
16697177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1674b9ad928SBarry Smith @*/
1687087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1694b9ad928SBarry Smith {
170dfbe8321SBarry Smith   PetscErrorCode ierr;
171f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
172f3fbd535SBarry Smith   MPI_Comm       comm = ((PetscObject)pc)->comm;
1733aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
174f3fbd535SBarry Smith   PetscInt       i;
175f3fbd535SBarry Smith   PetscMPIInt    size;
176f3fbd535SBarry Smith   const char     *prefix;
177f3fbd535SBarry Smith   PC             ipc;
1783aeaf226SBarry Smith   PetscInt       n;
1794b9ad928SBarry Smith 
1804b9ad928SBarry Smith   PetscFunctionBegin;
1810700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
182548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
183548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1843aeaf226SBarry Smith   if (mglevels) {
1853aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
1863aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
1873aeaf226SBarry Smith     n = mglevels[0]->levels;
1883aeaf226SBarry Smith     for (i=0; i<n; i++) {
1893aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1903aeaf226SBarry Smith 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
1913aeaf226SBarry Smith       }
1923aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
1933aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
1943aeaf226SBarry Smith     }
1953aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
1963aeaf226SBarry Smith   }
197f3fbd535SBarry Smith 
198f3fbd535SBarry Smith   mg->nlevels      = levels;
199f3fbd535SBarry Smith 
200f3fbd535SBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr);
201f3fbd535SBarry Smith   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
202f3fbd535SBarry Smith 
203f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
204f3fbd535SBarry Smith 
205a9db87a2SMatthew G Knepley   mg->stageApply = 0;
206f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
207f3fbd535SBarry Smith     ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr);
208f3fbd535SBarry Smith     mglevels[i]->level           = i;
209f3fbd535SBarry Smith     mglevels[i]->levels          = levels;
210f3fbd535SBarry Smith     mglevels[i]->cycles          = PC_MG_CYCLE_V;
211336babb1SJed Brown     mg->default_smoothu = 2;
212336babb1SJed Brown     mg->default_smoothd = 2;
21363e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
21463e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
21563e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
21663e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
217f3fbd535SBarry Smith 
218f3fbd535SBarry Smith     if (comms) comm = comms[i];
219f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
220336babb1SJed Brown     ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
221336babb1SJed Brown     ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
222336babb1SJed Brown     ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
223336babb1SJed Brown     ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
224336babb1SJed Brown     ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
225f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
226336babb1SJed Brown     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i?mg->default_smoothd:1);CHKERRQ(ierr);
227f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
228f3fbd535SBarry Smith 
229f3fbd535SBarry Smith     /* do special stuff for coarse grid */
230f3fbd535SBarry Smith     if (!i && levels > 1) {
231f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
232f3fbd535SBarry Smith 
233f3fbd535SBarry Smith       /* coarse solve is (redundant) LU by default */
234f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
235f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
236f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
237f3fbd535SBarry Smith       if (size > 1) {
238f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
239f3fbd535SBarry Smith       } else {
240f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
241f3fbd535SBarry Smith       }
242f3fbd535SBarry Smith 
243f3fbd535SBarry Smith     } else {
244f3fbd535SBarry Smith       char tprefix[128];
245f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
246f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
247f3fbd535SBarry Smith     }
248f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr);
249f3fbd535SBarry Smith     mglevels[i]->smoothu    = mglevels[i]->smoothd;
25031567311SBarry Smith     mg->rtol                = 0.0;
25131567311SBarry Smith     mg->abstol              = 0.0;
25231567311SBarry Smith     mg->dtol                = 0.0;
25331567311SBarry Smith     mg->ttol                = 0.0;
25431567311SBarry Smith     mg->cyclesperpcapply    = 1;
255f3fbd535SBarry Smith   }
25631567311SBarry Smith   mg->am          = PC_MG_MULTIPLICATIVE;
257f3fbd535SBarry Smith   mg->levels      = mglevels;
2584b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
2594b9ad928SBarry Smith   PetscFunctionReturn(0);
2604b9ad928SBarry Smith }
2614b9ad928SBarry Smith 
262c07bf074SBarry Smith 
263c07bf074SBarry Smith #undef __FUNCT__
264c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
265c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
266c07bf074SBarry Smith {
267c07bf074SBarry Smith   PetscErrorCode ierr;
268a06653b4SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
269a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
270a06653b4SBarry Smith   PetscInt       i,n;
271c07bf074SBarry Smith 
272c07bf074SBarry Smith   PetscFunctionBegin;
273a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
274a06653b4SBarry Smith   if (mglevels) {
275a06653b4SBarry Smith     n = mglevels[0]->levels;
276a06653b4SBarry Smith     for (i=0; i<n; i++) {
277a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2786bf464f9SBarry Smith 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
279a06653b4SBarry Smith       }
2806bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
281a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
282a06653b4SBarry Smith     }
283a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
284a06653b4SBarry Smith   }
285c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
286f3fbd535SBarry Smith   PetscFunctionReturn(0);
287f3fbd535SBarry Smith }
288f3fbd535SBarry Smith 
289f3fbd535SBarry Smith 
290f3fbd535SBarry Smith 
29109573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
29209573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
29309573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
294f3fbd535SBarry Smith 
295f3fbd535SBarry Smith /*
296f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
297f3fbd535SBarry Smith              or full cycle of multigrid.
298f3fbd535SBarry Smith 
299f3fbd535SBarry Smith   Note:
300f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
301f3fbd535SBarry Smith */
302f3fbd535SBarry Smith #undef __FUNCT__
303f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
304f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
305f3fbd535SBarry Smith {
306f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
307f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
308f3fbd535SBarry Smith   PetscErrorCode ierr;
309f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
310f3fbd535SBarry Smith 
311f3fbd535SBarry Smith   PetscFunctionBegin;
312a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
313e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
314298cc208SBarry Smith   for (i=0; i<levels; i++) {
315e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
316e1d8e5deSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
317298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
318e1d8e5deSBarry Smith     }
319e1d8e5deSBarry Smith   }
320e1d8e5deSBarry Smith 
321f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
322f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
32331567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
324f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
32531567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
326f3fbd535SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
327f3fbd535SBarry Smith     }
328f3fbd535SBarry Smith   }
32931567311SBarry Smith   else if (mg->am == PC_MG_ADDITIVE) {
33031567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
331f3fbd535SBarry Smith   }
33231567311SBarry Smith   else if (mg->am == PC_MG_KASKADE) {
33331567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
334f3fbd535SBarry Smith   }
335f3fbd535SBarry Smith   else {
336f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
337f3fbd535SBarry Smith   }
338a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
339f3fbd535SBarry Smith   PetscFunctionReturn(0);
340f3fbd535SBarry Smith }
341f3fbd535SBarry Smith 
342f3fbd535SBarry Smith 
343f3fbd535SBarry Smith #undef __FUNCT__
344f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
345f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc)
346f3fbd535SBarry Smith {
347f3fbd535SBarry Smith   PetscErrorCode ierr;
348f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
349c91913d3SJed Brown   PetscBool      flg,set;
350f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
351f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
352f3fbd535SBarry Smith   PCMGType       mgtype;
353f3fbd535SBarry Smith   PCMGCycleType  mgctype;
354f3fbd535SBarry Smith 
355f3fbd535SBarry Smith   PetscFunctionBegin;
356f3fbd535SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
3573aeaf226SBarry Smith     if (!mg->levels) {
358f3fbd535SBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
359eb3f98d2SBarry Smith       if (!flg && pc->dm) {
360eb3f98d2SBarry Smith         ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
361eb3f98d2SBarry Smith         levels++;
3623aeaf226SBarry Smith         mg->usedmfornumberoflevels = PETSC_TRUE;
363eb3f98d2SBarry Smith       }
364f3fbd535SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
365f3fbd535SBarry Smith     }
3663aeaf226SBarry Smith     mglevels = mg->levels;
3673aeaf226SBarry Smith 
368f3fbd535SBarry Smith     mgctype = (PCMGCycleType) mglevels[0]->cycles;
369f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
370f3fbd535SBarry Smith     if (flg) {
371f3fbd535SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
372f3fbd535SBarry Smith     };
373f3fbd535SBarry Smith     flg  = PETSC_FALSE;
374c91913d3SJed Brown     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr);
375c91913d3SJed Brown     if (set) {
376c91913d3SJed Brown       ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr);
377f3fbd535SBarry Smith     }
378f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
379f3fbd535SBarry Smith     if (flg) {
380f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
381f3fbd535SBarry Smith     }
382f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
383f3fbd535SBarry Smith     if (flg) {
384f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
385f3fbd535SBarry Smith     }
38631567311SBarry Smith     mgtype = mg->am;
387f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
388f3fbd535SBarry Smith     if (flg) {
389f3fbd535SBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
390f3fbd535SBarry Smith     }
39131567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
39231567311SBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
393f3fbd535SBarry Smith       if (flg) {
394f3fbd535SBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
395f3fbd535SBarry Smith       }
396f3fbd535SBarry Smith     }
397f3fbd535SBarry Smith     flg  = PETSC_FALSE;
398acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
399f3fbd535SBarry Smith     if (flg) {
400f3fbd535SBarry Smith       PetscInt i;
401f3fbd535SBarry Smith       char     eventname[128];
40263e6d426SJed Brown       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
403f3fbd535SBarry Smith       levels = mglevels[0]->levels;
404f3fbd535SBarry Smith       for (i=0; i<levels; i++) {
405f3fbd535SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
40663e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
407f3fbd535SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
40863e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
409f3fbd535SBarry Smith         if (i) {
410f3fbd535SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
41163e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
412f3fbd535SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
41363e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
414f3fbd535SBarry Smith         }
415f3fbd535SBarry Smith       }
416eec35531SMatthew G Knepley 
417ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
418eec35531SMatthew G Knepley       {
419eec35531SMatthew G Knepley         const char   *sname = "MG Apply";
420eec35531SMatthew G Knepley         PetscStageLog stageLog;
421eec35531SMatthew G Knepley         PetscInt      st;
422eec35531SMatthew G Knepley 
423eec35531SMatthew G Knepley         PetscFunctionBegin;
424eec35531SMatthew G Knepley         ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
425eec35531SMatthew G Knepley         for (st = 0; st < stageLog->numStages; ++st) {
426eec35531SMatthew G Knepley           PetscBool same;
427eec35531SMatthew G Knepley 
428eec35531SMatthew G Knepley           ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
429eec35531SMatthew G Knepley           if (same) {mg->stageApply = st;}
430eec35531SMatthew G Knepley         }
431eec35531SMatthew G Knepley         if (!mg->stageApply) {
432eec35531SMatthew G Knepley           ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
433eec35531SMatthew G Knepley         }
434eec35531SMatthew G Knepley       }
435ec60ca73SSatish Balay #endif
436f3fbd535SBarry Smith     }
437f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
438f3fbd535SBarry Smith   PetscFunctionReturn(0);
439f3fbd535SBarry Smith }
440f3fbd535SBarry Smith 
4416a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4426a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
443f3fbd535SBarry Smith 
444f3fbd535SBarry Smith #undef __FUNCT__
445f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
446f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
447f3fbd535SBarry Smith {
448f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
449f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
450f3fbd535SBarry Smith   PetscErrorCode ierr;
451f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
452219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
453f3fbd535SBarry Smith 
454f3fbd535SBarry Smith   PetscFunctionBegin;
455251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4565b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
457219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
458f3fbd535SBarry Smith   if (iascii) {
45931567311SBarry 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);
46031567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
46131567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
462f3fbd535SBarry Smith     }
463218a07d4SBarry Smith     if (mg->galerkin) {
464f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4654f66f45eSBarry Smith     } else {
4664f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
467f3fbd535SBarry Smith     }
468f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
469f3fbd535SBarry Smith       if (!i) {
47089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
471f3fbd535SBarry Smith       } else {
47289cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
473f3fbd535SBarry Smith       }
474f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
475f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
476f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
477f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
478f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
479f3fbd535SBarry Smith       } else if (i){
48089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
481f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
482f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
483f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
484f3fbd535SBarry Smith       }
485f3fbd535SBarry Smith     }
4865b0b0462SBarry Smith   } else if (isbinary) {
4875b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
4885b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
4895b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
4905b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
4915b0b0462SBarry Smith       }
4925b0b0462SBarry Smith     }
493219b1045SBarry Smith   } else if (isdraw) {
494219b1045SBarry Smith     PetscDraw draw;
495219b1045SBarry Smith     PetscReal x,y,bottom,th;
496219b1045SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
497219b1045SBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
498219b1045SBarry Smith     ierr = PetscDrawStringGetSize(draw,PETSC_NULL,&th);CHKERRQ(ierr);
499219b1045SBarry Smith     bottom = y - th;
500219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
501219b1045SBarry Smith       ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
502219b1045SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
503219b1045SBarry Smith       ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
504*72a02f06SBarry Smith       /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
505219b1045SBarry Smith       bottom -= 5*th;
506219b1045SBarry Smith     }
5075b0b0462SBarry Smith   }
508f3fbd535SBarry Smith   PetscFunctionReturn(0);
509f3fbd535SBarry Smith }
510f3fbd535SBarry Smith 
511b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>
512b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
513cab2e9ccSBarry Smith 
514f3fbd535SBarry Smith /*
515f3fbd535SBarry Smith     Calls setup for the KSP on each level
516f3fbd535SBarry Smith */
517f3fbd535SBarry Smith #undef __FUNCT__
518f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
519f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
520f3fbd535SBarry Smith {
521f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
522f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
523f3fbd535SBarry Smith   PetscErrorCode          ierr;
524f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
52598e04984SBarry Smith   PC                      cpc;
526649052a6SBarry Smith   PetscBool               preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset;
527f3fbd535SBarry Smith   Mat                     dA,dB;
528f3fbd535SBarry Smith   MatStructure            uflag;
529f3fbd535SBarry Smith   Vec                     tvec;
530218a07d4SBarry Smith   DM                      *dms;
531649052a6SBarry Smith   PetscViewer             viewer = 0;
532f3fbd535SBarry Smith 
533f3fbd535SBarry Smith   PetscFunctionBegin;
53401bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5353aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5363aeaf226SBarry Smith     PetscInt levels;
5373aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5383aeaf226SBarry Smith     levels++;
5393aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5403aeaf226SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
5413aeaf226SBarry Smith       n    = levels;
5423aeaf226SBarry Smith       ierr = PCSetFromOptions(pc);CHKERRQ(ierr);  /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5433aeaf226SBarry Smith       mglevels =  mg->levels;
5443aeaf226SBarry Smith     }
5453aeaf226SBarry Smith   }
54698e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5473aeaf226SBarry Smith 
548f3fbd535SBarry Smith 
549f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
550f3fbd535SBarry Smith   /* so use those from global PC */
551f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
552f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
55398e04984SBarry Smith   if (opsset) {
55498e04984SBarry Smith     Mat mmat;
55598e04984SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,PETSC_NULL,&mmat,PETSC_NULL);CHKERRQ(ierr);
55698e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
55798e04984SBarry Smith   }
55898e04984SBarry Smith   if (!opsset) {
559f3fbd535SBarry 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);
560f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
561f3fbd535SBarry Smith   }
562f3fbd535SBarry Smith 
56301bc414fSDmitry Karpeev   /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */
564ce4cda84SJed Brown   if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
5652d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
5662e499ae9SBarry Smith     Mat p;
56773250ac0SBarry Smith     Vec rscale;
568218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
569218a07d4SBarry Smith     dms[n-1] = pc->dm;
570218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
571942e3340SBarry Smith       DMKSP kdm;
5722ee06e3bSJed Brown       ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);
5733c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
5743c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
575942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
576d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
577d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
5785358d0d4SBarry Smith       kdm->ops->computerhs = PETSC_NULL;
579fe86f630SJed Brown       kdm->rhsctx          = PETSC_NULL;
58024c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
581e727c939SJed Brown 	ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
5822d2b81a6SBarry Smith 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
58300ac8be1SBarry Smith 	if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
58473250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
5856bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
586218a07d4SBarry Smith       }
58724c3aa18SBarry Smith     }
5882d2b81a6SBarry Smith 
589218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
5906bf464f9SBarry Smith       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
591218a07d4SBarry Smith     }
5922d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
593ce4cda84SJed Brown   }
594cab2e9ccSBarry Smith 
595ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
596ce4cda84SJed Brown     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
597cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
598cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
599218a07d4SBarry Smith   }
600218a07d4SBarry Smith 
601c91913d3SJed Brown   if (mg->galerkin == 1) {
602f3fbd535SBarry Smith     Mat B;
603f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
604f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
605f3fbd535SBarry Smith     if (!pc->setupcalled) {
606f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
607f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
608f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
609f3fbd535SBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
610f3fbd535SBarry Smith         dB   = B;
611f3fbd535SBarry Smith       }
612cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
613f3fbd535SBarry Smith     } else {
614f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
615f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
616f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
617f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
618f3fbd535SBarry Smith         dB   = B;
619f3fbd535SBarry Smith       }
620f3fbd535SBarry Smith     }
621ce4cda84SJed Brown   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
622cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
623cab2e9ccSBarry Smith     for (i=n-2;i>-1; i--) {
624c88c5224SJed Brown       Mat R;
625c88c5224SJed Brown       Vec rscale;
626cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
627cab2e9ccSBarry Smith         Vec *vecs;
628cab2e9ccSBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr);
629cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
630cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
631cab2e9ccSBarry Smith       }
632c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
633c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
634c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
635c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
636cab2e9ccSBarry Smith     }
637f3fbd535SBarry Smith   }
638ccceb50cSJed Brown   if (!mg->galerkin && pc->dm) {
639caa4e7f2SJed Brown     for (i=n-2;i>=0; i--) {
640ccceb50cSJed Brown       DM dmfine,dmcoarse;
641caa4e7f2SJed Brown       Mat Restrict,Inject;
642caa4e7f2SJed Brown       Vec rscale;
643ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
644ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
645caa4e7f2SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
646caa4e7f2SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
647caa4e7f2SJed Brown       Inject = PETSC_NULL;      /* Callback should create it if it needs Injection */
648caa4e7f2SJed Brown       ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
649caa4e7f2SJed Brown     }
650caa4e7f2SJed Brown   }
651f3fbd535SBarry Smith 
652f3fbd535SBarry Smith   if (!pc->setupcalled) {
653f3fbd535SBarry Smith     for (i=0; i<n; i++) {
654f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
655f3fbd535SBarry Smith     }
656f3fbd535SBarry Smith     for (i=1; i<n; i++) {
657f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
658f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
659f3fbd535SBarry Smith       }
660f3fbd535SBarry Smith     }
661f3fbd535SBarry Smith     for (i=1; i<n; i++) {
662c88c5224SJed Brown       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
663c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
664f3fbd535SBarry Smith     }
665f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
666f3fbd535SBarry Smith       if (!mglevels[i]->b) {
667f3fbd535SBarry Smith         Vec *vec;
668f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
669f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
6706bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
671f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
672f3fbd535SBarry Smith       }
673f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
674f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
675f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
6766bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
677f3fbd535SBarry Smith       }
678f3fbd535SBarry Smith       if (!mglevels[i]->x) {
679f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
680f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
6816bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
682f3fbd535SBarry Smith       }
683f3fbd535SBarry Smith     }
684f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
685f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
686f3fbd535SBarry Smith       Vec *vec;
687f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
688f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
6896bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
690f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
691f3fbd535SBarry Smith     }
692f3fbd535SBarry Smith   }
693f3fbd535SBarry Smith 
69498e04984SBarry Smith   if (pc->dm) {
69598e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
69698e04984SBarry Smith     for (i=0; i<n-1; i++){
69798e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
69898e04984SBarry Smith     }
69998e04984SBarry Smith   }
700f3fbd535SBarry Smith 
701f3fbd535SBarry Smith   for (i=1; i<n; i++) {
702f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
703f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
704f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
705f3fbd535SBarry Smith     }
70663e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
707f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
70863e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
709d42688cbSBarry Smith     if (!mglevels[i]->residual) {
710d42688cbSBarry Smith       Mat mat;
711d42688cbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
712d42688cbSBarry Smith       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
713d42688cbSBarry Smith     }
714f3fbd535SBarry Smith   }
715f3fbd535SBarry Smith   for (i=1; i<n; i++) {
716f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
717f3fbd535SBarry Smith       Mat          downmat,downpmat;
718f3fbd535SBarry Smith       MatStructure matflag;
719f3fbd535SBarry Smith 
720f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
721f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
722f3fbd535SBarry Smith       if (!opsset) {
723f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
724f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
725f3fbd535SBarry Smith       }
726f3fbd535SBarry Smith 
727f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
72863e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
729f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
73063e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
731f3fbd535SBarry Smith     }
732f3fbd535SBarry Smith   }
733f3fbd535SBarry Smith 
734f3fbd535SBarry Smith   /*
735f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
736f3fbd535SBarry Smith   */
737251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
738f3fbd535SBarry Smith   if (preonly) {
739251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
740251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
741251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
742251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
743fd1303e9SJungho Lee     if (!lu && !redundant && !cholesky && !svd) {
744f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
745f3fbd535SBarry Smith     }
746f3fbd535SBarry Smith   }
747f3fbd535SBarry Smith 
748f3fbd535SBarry Smith   if (!pc->setupcalled) {
749f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
750f3fbd535SBarry Smith   }
751f3fbd535SBarry Smith 
75263e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
753f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
75463e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
755f3fbd535SBarry Smith 
756f3fbd535SBarry Smith   /*
757f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
758e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
759f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
760f3fbd535SBarry Smith 
761f3fbd535SBarry Smith    Only support one or the other at the same time.
762f3fbd535SBarry Smith   */
763f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
764acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
765f3fbd535SBarry Smith   if (dump) {
766f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
767f3fbd535SBarry Smith   }
768f3fbd535SBarry Smith   dump = PETSC_FALSE;
769f3fbd535SBarry Smith #endif
770acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
771f3fbd535SBarry Smith   if (dump) {
77232c0f0efSBarry Smith 
773f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
774f3fbd535SBarry Smith   }
775f3fbd535SBarry Smith 
776f3fbd535SBarry Smith   if (viewer) {
777f3fbd535SBarry Smith     for (i=1; i<n; i++) {
778f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
779f3fbd535SBarry Smith     }
780f3fbd535SBarry Smith     for (i=0; i<n; i++) {
781f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
782f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
783f3fbd535SBarry Smith     }
784f3fbd535SBarry Smith   }
785f3fbd535SBarry Smith   PetscFunctionReturn(0);
786f3fbd535SBarry Smith }
787f3fbd535SBarry Smith 
788f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
789f3fbd535SBarry Smith 
790f3fbd535SBarry Smith #undef __FUNCT__
7919dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
7924b9ad928SBarry Smith /*@
79397177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
7944b9ad928SBarry Smith 
7954b9ad928SBarry Smith    Not Collective
7964b9ad928SBarry Smith 
7974b9ad928SBarry Smith    Input Parameter:
7984b9ad928SBarry Smith .  pc - the preconditioner context
7994b9ad928SBarry Smith 
8004b9ad928SBarry Smith    Output parameter:
8014b9ad928SBarry Smith .  levels - the number of levels
8024b9ad928SBarry Smith 
8034b9ad928SBarry Smith    Level: advanced
8044b9ad928SBarry Smith 
8054b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8064b9ad928SBarry Smith 
80797177400SBarry Smith .seealso: PCMGSetLevels()
8084b9ad928SBarry Smith @*/
8097087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8104b9ad928SBarry Smith {
811f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8124b9ad928SBarry Smith 
8134b9ad928SBarry Smith   PetscFunctionBegin;
8140700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8154482741eSBarry Smith   PetscValidIntPointer(levels,2);
816f3fbd535SBarry Smith   *levels = mg->nlevels;
8174b9ad928SBarry Smith   PetscFunctionReturn(0);
8184b9ad928SBarry Smith }
8194b9ad928SBarry Smith 
8204b9ad928SBarry Smith #undef __FUNCT__
8219dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
8224b9ad928SBarry Smith /*@
82397177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8244b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8254b9ad928SBarry Smith 
826ad4df100SBarry Smith    Logically Collective on PC
8274b9ad928SBarry Smith 
8284b9ad928SBarry Smith    Input Parameters:
8294b9ad928SBarry Smith +  pc - the preconditioner context
8309dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
8319dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
8324b9ad928SBarry Smith 
8334b9ad928SBarry Smith    Options Database Key:
8344b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
8354b9ad928SBarry Smith    additive, full, kaskade
8364b9ad928SBarry Smith 
8374b9ad928SBarry Smith    Level: advanced
8384b9ad928SBarry Smith 
8394b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
8404b9ad928SBarry Smith 
84197177400SBarry Smith .seealso: PCMGSetLevels()
8424b9ad928SBarry Smith @*/
8437087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
8444b9ad928SBarry Smith {
845f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
8464b9ad928SBarry Smith 
8474b9ad928SBarry Smith   PetscFunctionBegin;
8480700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
849c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
85031567311SBarry Smith   mg->am = form;
8519dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
8524b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
8534b9ad928SBarry Smith   PetscFunctionReturn(0);
8544b9ad928SBarry Smith }
8554b9ad928SBarry Smith 
8564b9ad928SBarry Smith #undef __FUNCT__
8570d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
8584b9ad928SBarry Smith /*@
8590d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
8604b9ad928SBarry Smith    complicated cycling.
8614b9ad928SBarry Smith 
862ad4df100SBarry Smith    Logically Collective on PC
8634b9ad928SBarry Smith 
8644b9ad928SBarry Smith    Input Parameters:
865c2be2410SBarry Smith +  pc - the multigrid context
8660d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
8674b9ad928SBarry Smith 
8684b9ad928SBarry Smith    Options Database Key:
8690d353602SBarry Smith $  -pc_mg_cycle_type v or w
8704b9ad928SBarry Smith 
8714b9ad928SBarry Smith    Level: advanced
8724b9ad928SBarry Smith 
8734b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8744b9ad928SBarry Smith 
8750d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
8764b9ad928SBarry Smith @*/
8777087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
8784b9ad928SBarry Smith {
879f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
880f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
88179416396SBarry Smith   PetscInt     i,levels;
8824b9ad928SBarry Smith 
8834b9ad928SBarry Smith   PetscFunctionBegin;
8840700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
885e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
886c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
887f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8884b9ad928SBarry Smith 
8894b9ad928SBarry Smith   for (i=0; i<levels; i++) {
890f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
8914b9ad928SBarry Smith   }
8924b9ad928SBarry Smith   PetscFunctionReturn(0);
8934b9ad928SBarry Smith }
8944b9ad928SBarry Smith 
8954b9ad928SBarry Smith #undef __FUNCT__
8968cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
8978cc2d5dfSBarry Smith /*@
8988cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
8998cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9008cc2d5dfSBarry Smith 
901ad4df100SBarry Smith    Logically Collective on PC
9028cc2d5dfSBarry Smith 
9038cc2d5dfSBarry Smith    Input Parameters:
9048cc2d5dfSBarry Smith +  pc - the multigrid context
9058cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
9068cc2d5dfSBarry Smith 
9078cc2d5dfSBarry Smith    Options Database Key:
9088cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
9098cc2d5dfSBarry Smith 
9108cc2d5dfSBarry Smith    Level: advanced
9118cc2d5dfSBarry Smith 
9128cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
9138cc2d5dfSBarry Smith 
9148cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9158cc2d5dfSBarry Smith 
9168cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
9178cc2d5dfSBarry Smith @*/
9187087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
9198cc2d5dfSBarry Smith {
920f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
921f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
9228cc2d5dfSBarry Smith   PetscInt     i,levels;
9238cc2d5dfSBarry Smith 
9248cc2d5dfSBarry Smith   PetscFunctionBegin;
9250700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
926e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
927c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
928f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9298cc2d5dfSBarry Smith 
9308cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
93131567311SBarry Smith     mg->cyclesperpcapply  = n;
9328cc2d5dfSBarry Smith   }
9338cc2d5dfSBarry Smith   PetscFunctionReturn(0);
9348cc2d5dfSBarry Smith }
9358cc2d5dfSBarry Smith 
9368cc2d5dfSBarry Smith #undef __FUNCT__
9379dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
938c2be2410SBarry Smith /*@
93997177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
940c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
941c2be2410SBarry Smith 
942ad4df100SBarry Smith    Logically Collective on PC
943c2be2410SBarry Smith 
944c2be2410SBarry Smith    Input Parameters:
945c91913d3SJed Brown +  pc - the multigrid context
946c91913d3SJed Brown -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
947c2be2410SBarry Smith 
948c2be2410SBarry Smith    Options Database Key:
949c2be2410SBarry Smith $  -pc_mg_galerkin
950c2be2410SBarry Smith 
951c2be2410SBarry Smith    Level: intermediate
952c2be2410SBarry Smith 
953c2be2410SBarry Smith .keywords: MG, set, Galerkin
954c2be2410SBarry Smith 
9553fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
9563fc8bf9cSBarry Smith 
957c2be2410SBarry Smith @*/
958c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
959c2be2410SBarry Smith {
960f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
961c2be2410SBarry Smith 
962c2be2410SBarry Smith   PetscFunctionBegin;
9630700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
964789726d7SBarry Smith   mg->galerkin = use ? 1 : 0;
965c2be2410SBarry Smith   PetscFunctionReturn(0);
966c2be2410SBarry Smith }
967c2be2410SBarry Smith 
968c2be2410SBarry Smith #undef __FUNCT__
9693fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
9703fc8bf9cSBarry Smith /*@
9713fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
9723fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
9733fc8bf9cSBarry Smith 
9743fc8bf9cSBarry Smith    Not Collective
9753fc8bf9cSBarry Smith 
9763fc8bf9cSBarry Smith    Input Parameter:
9773fc8bf9cSBarry Smith .  pc - the multigrid context
9783fc8bf9cSBarry Smith 
9793fc8bf9cSBarry Smith    Output Parameter:
9803fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
9813fc8bf9cSBarry Smith 
9823fc8bf9cSBarry Smith    Options Database Key:
9833fc8bf9cSBarry Smith $  -pc_mg_galerkin
9843fc8bf9cSBarry Smith 
9853fc8bf9cSBarry Smith    Level: intermediate
9863fc8bf9cSBarry Smith 
9873fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
9883fc8bf9cSBarry Smith 
9893fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
9903fc8bf9cSBarry Smith 
9913fc8bf9cSBarry Smith @*/
9927087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
9933fc8bf9cSBarry Smith {
994f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
9953fc8bf9cSBarry Smith 
9963fc8bf9cSBarry Smith   PetscFunctionBegin;
9970700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
99813fdb3acSJose Roman   *galerkin = (PetscBool)mg->galerkin;
9993fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10003fc8bf9cSBarry Smith }
10013fc8bf9cSBarry Smith 
10023fc8bf9cSBarry Smith #undef __FUNCT__
10039dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
10044b9ad928SBarry Smith /*@
100597177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
100697177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
10074b9ad928SBarry Smith    pre-smoothing steps on different levels.
10084b9ad928SBarry Smith 
1009ad4df100SBarry Smith    Logically Collective on PC
10104b9ad928SBarry Smith 
10114b9ad928SBarry Smith    Input Parameters:
10124b9ad928SBarry Smith +  mg - the multigrid context
10134b9ad928SBarry Smith -  n - the number of smoothing steps
10144b9ad928SBarry Smith 
10154b9ad928SBarry Smith    Options Database Key:
10164b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
10174b9ad928SBarry Smith 
10184b9ad928SBarry Smith    Level: advanced
10194b9ad928SBarry Smith 
10204b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
10214b9ad928SBarry Smith 
102297177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
10234b9ad928SBarry Smith @*/
10247087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
10254b9ad928SBarry Smith {
1026f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1027f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10286849ba73SBarry Smith   PetscErrorCode ierr;
102979416396SBarry Smith   PetscInt       i,levels;
10304b9ad928SBarry Smith 
10314b9ad928SBarry Smith   PetscFunctionBegin;
10320700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1033e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1034c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1035f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10364b9ad928SBarry Smith 
1037b05257ddSBarry Smith   for (i=1; i<levels; i++) {
10384b9ad928SBarry Smith     /* make sure smoother up and down are different */
103997177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1040f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
104131567311SBarry Smith     mg->default_smoothd = n;
10424b9ad928SBarry Smith   }
10434b9ad928SBarry Smith   PetscFunctionReturn(0);
10444b9ad928SBarry Smith }
10454b9ad928SBarry Smith 
10464b9ad928SBarry Smith #undef __FUNCT__
10479dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
10484b9ad928SBarry Smith /*@
104997177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
105097177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
10514b9ad928SBarry Smith    post-smoothing steps on different levels.
10524b9ad928SBarry Smith 
1053ad4df100SBarry Smith    Logically Collective on PC
10544b9ad928SBarry Smith 
10554b9ad928SBarry Smith    Input Parameters:
10564b9ad928SBarry Smith +  mg - the multigrid context
10574b9ad928SBarry Smith -  n - the number of smoothing steps
10584b9ad928SBarry Smith 
10594b9ad928SBarry Smith    Options Database Key:
10604b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
10614b9ad928SBarry Smith 
10624b9ad928SBarry Smith    Level: advanced
10634b9ad928SBarry Smith 
10644b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
1065a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
10664b9ad928SBarry Smith 
10674b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
10684b9ad928SBarry Smith 
106997177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
10704b9ad928SBarry Smith @*/
10717087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
10724b9ad928SBarry Smith {
1073f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1074f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10756849ba73SBarry Smith   PetscErrorCode ierr;
107679416396SBarry Smith   PetscInt       i,levels;
10774b9ad928SBarry Smith 
10784b9ad928SBarry Smith   PetscFunctionBegin;
10790700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1080e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1081c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1082f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10834b9ad928SBarry Smith 
10844b9ad928SBarry Smith   for (i=1; i<levels; i++) {
10854b9ad928SBarry Smith     /* make sure smoother up and down are different */
108697177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1087f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
108831567311SBarry Smith     mg->default_smoothu = n;
10894b9ad928SBarry Smith   }
10904b9ad928SBarry Smith   PetscFunctionReturn(0);
10914b9ad928SBarry Smith }
10924b9ad928SBarry Smith 
10934b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
10944b9ad928SBarry Smith 
10953b09bd56SBarry Smith /*MC
1096ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
10973b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
10983b09bd56SBarry Smith 
10993b09bd56SBarry Smith    Options Database Keys:
11003b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
11010d353602SBarry Smith .  -pc_mg_cycles v or w
110279416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
11033b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
11048c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
11053b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
11063b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
11078cf6037eSBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
11088cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
11098cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1110e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
11118cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
11128cf6037eSBarry Smith                         to the binary output file called binaryoutput
11133b09bd56SBarry Smith 
111424c3aa18SBarry 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
11153b09bd56SBarry Smith 
11168cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
11178cf6037eSBarry Smith 
11183b09bd56SBarry Smith    Level: intermediate
11193b09bd56SBarry Smith 
11208f87f92bSBarry Smith    Concepts: multigrid/multilevel
11213b09bd56SBarry Smith 
11228cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
11230d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
112497177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
112597177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
11260d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11273b09bd56SBarry Smith M*/
11283b09bd56SBarry Smith 
11294b9ad928SBarry Smith EXTERN_C_BEGIN
11304b9ad928SBarry Smith #undef __FUNCT__
11314b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
11327087cfbeSBarry Smith PetscErrorCode  PCCreate_MG(PC pc)
11334b9ad928SBarry Smith {
1134f3fbd535SBarry Smith   PC_MG          *mg;
1135f3fbd535SBarry Smith   PetscErrorCode ierr;
1136f3fbd535SBarry Smith 
11374b9ad928SBarry Smith   PetscFunctionBegin;
1138f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1139f3fbd535SBarry Smith   pc->data    = (void*)mg;
1140f3fbd535SBarry Smith   mg->nlevels = -1;
1141f3fbd535SBarry Smith 
11424b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
11434b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1144a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
11454b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
11464b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
11474b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
11484b9ad928SBarry Smith   PetscFunctionReturn(0);
11494b9ad928SBarry Smith }
11504b9ad928SBarry Smith EXTERN_C_END
1151