xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 00ac8be1419b46f104fcdb8335f77affb930ddf0)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscpcmg.h" I*/
64b9ad928SBarry Smith 
74b9ad928SBarry Smith 
84b9ad928SBarry Smith #undef __FUNCT__
99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private"
1031567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
114b9ad928SBarry Smith {
1231567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1331567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
146849ba73SBarry Smith   PetscErrorCode ierr;
1531567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
184b9ad928SBarry Smith 
1963e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2031567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
2163e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2231567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2363e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2431567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2563e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
264b9ad928SBarry Smith 
274b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
2831567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
294b9ad928SBarry Smith       PetscReal rnorm;
3031567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
314b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3270441072SBarry Smith         if (rnorm < mg->abstol) {
334d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
341e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
354b9ad928SBarry Smith         } else {
364d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
371e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr);
384b9ad928SBarry Smith         }
394b9ad928SBarry Smith         PetscFunctionReturn(0);
404b9ad928SBarry Smith       }
414b9ad928SBarry Smith     }
424b9ad928SBarry Smith 
4331567311SBarry Smith     mgc = *(mglevelsin - 1);
4463e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4531567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
4663e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
47efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
484b9ad928SBarry Smith     while (cycles--) {
4931567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
504b9ad928SBarry Smith     }
5163e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5231567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5363e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5463e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5531567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
5663e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith   }
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
614b9ad928SBarry Smith #undef __FUNCT__
624b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
63ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool  zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
644b9ad928SBarry Smith {
65f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
66f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
67dfbe8321SBarry Smith   PetscErrorCode ierr;
68f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
694b9ad928SBarry Smith 
704b9ad928SBarry Smith   PetscFunctionBegin;
71f3fbd535SBarry Smith   mglevels[levels-1]->b    = b;
72f3fbd535SBarry Smith   mglevels[levels-1]->x    = x;
734b9ad928SBarry Smith 
7431567311SBarry Smith   mg->rtol = rtol;
7531567311SBarry Smith   mg->abstol = abstol;
7631567311SBarry Smith   mg->dtol = dtol;
774b9ad928SBarry Smith   if (rtol) {
784b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
794b9ad928SBarry Smith     PetscReal rnorm;
807319c654SBarry Smith     if (zeroguess) {
817319c654SBarry Smith       ierr               = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
827319c654SBarry Smith     } else {
83f3fbd535SBarry Smith       ierr               = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
844b9ad928SBarry Smith       ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
857319c654SBarry Smith     }
8631567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
8770441072SBarry Smith   } else if (abstol) {
8831567311SBarry Smith     mg->ttol = abstol;
894b9ad928SBarry Smith   } else {
9031567311SBarry Smith     mg->ttol = 0.0;
914b9ad928SBarry Smith   }
924b9ad928SBarry Smith 
934d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
944d0a8057SBarry Smith      stop prematurely do to small residual */
954d0a8057SBarry Smith   for (i=1; i<levels; i++) {
96f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
97f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
98f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
994b9ad928SBarry Smith     }
1004d0a8057SBarry Smith   }
1014d0a8057SBarry Smith 
1024d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1034d0a8057SBarry Smith   for (i=0; i<its; i++) {
104f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1054d0a8057SBarry Smith     if (*reason) break;
1064d0a8057SBarry Smith   }
1074d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1084d0a8057SBarry Smith   *outits = i;
1094b9ad928SBarry Smith   PetscFunctionReturn(0);
1104b9ad928SBarry Smith }
1114b9ad928SBarry Smith 
1124b9ad928SBarry Smith #undef __FUNCT__
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;
199218a07d4SBarry Smith   mg->galerkin     = PETSC_FALSE;
200218a07d4SBarry Smith   mg->galerkinused = PETSC_FALSE;
201f3fbd535SBarry Smith 
202f3fbd535SBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr);
203f3fbd535SBarry Smith   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
204f3fbd535SBarry Smith 
205f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
206f3fbd535SBarry Smith 
207f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
208f3fbd535SBarry Smith     ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr);
209f3fbd535SBarry Smith     mglevels[i]->level           = i;
210f3fbd535SBarry Smith     mglevels[i]->levels          = levels;
211f3fbd535SBarry Smith     mglevels[i]->cycles          = PC_MG_CYCLE_V;
21231567311SBarry Smith     mg->default_smoothu = 1;
21331567311SBarry Smith     mg->default_smoothd = 1;
21463e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
21563e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
21663e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
21763e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
218f3fbd535SBarry Smith 
219f3fbd535SBarry Smith     if (comms) comm = comms[i];
220f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
221f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
22231567311SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
223f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
224f3fbd535SBarry Smith 
225f3fbd535SBarry Smith     /* do special stuff for coarse grid */
226f3fbd535SBarry Smith     if (!i && levels > 1) {
227f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
228f3fbd535SBarry Smith 
229f3fbd535SBarry Smith       /* coarse solve is (redundant) LU by default */
230f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
231f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
232f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
233f3fbd535SBarry Smith       if (size > 1) {
234f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
235f3fbd535SBarry Smith       } else {
236f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
237f3fbd535SBarry Smith       }
238f3fbd535SBarry Smith 
239f3fbd535SBarry Smith     } else {
240f3fbd535SBarry Smith       char tprefix[128];
241f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
242f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
243f3fbd535SBarry Smith     }
244f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr);
245f3fbd535SBarry Smith     mglevels[i]->smoothu    = mglevels[i]->smoothd;
24631567311SBarry Smith     mg->rtol                = 0.0;
24731567311SBarry Smith     mg->abstol              = 0.0;
24831567311SBarry Smith     mg->dtol                = 0.0;
24931567311SBarry Smith     mg->ttol                = 0.0;
25031567311SBarry Smith     mg->cyclesperpcapply    = 1;
251f3fbd535SBarry Smith   }
25231567311SBarry Smith   mg->am          = PC_MG_MULTIPLICATIVE;
253f3fbd535SBarry Smith   mg->levels      = mglevels;
2544b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
2554b9ad928SBarry Smith   PetscFunctionReturn(0);
2564b9ad928SBarry Smith }
2574b9ad928SBarry Smith 
258c07bf074SBarry Smith 
259c07bf074SBarry Smith #undef __FUNCT__
260c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
261c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
262c07bf074SBarry Smith {
263c07bf074SBarry Smith   PetscErrorCode ierr;
264a06653b4SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
265a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
266a06653b4SBarry Smith   PetscInt       i,n;
267c07bf074SBarry Smith 
268c07bf074SBarry Smith   PetscFunctionBegin;
269a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
270a06653b4SBarry Smith   if (mglevels) {
271a06653b4SBarry Smith     n = mglevels[0]->levels;
272a06653b4SBarry Smith     for (i=0; i<n; i++) {
273a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2746bf464f9SBarry Smith 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
275a06653b4SBarry Smith       }
2766bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
277a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
278a06653b4SBarry Smith     }
279a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
280a06653b4SBarry Smith   }
281c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
282f3fbd535SBarry Smith   PetscFunctionReturn(0);
283f3fbd535SBarry Smith }
284f3fbd535SBarry Smith 
285f3fbd535SBarry Smith 
286f3fbd535SBarry Smith 
28709573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
28809573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
28909573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
290f3fbd535SBarry Smith 
291f3fbd535SBarry Smith /*
292f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
293f3fbd535SBarry Smith              or full cycle of multigrid.
294f3fbd535SBarry Smith 
295f3fbd535SBarry Smith   Note:
296f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
297f3fbd535SBarry Smith */
298f3fbd535SBarry Smith #undef __FUNCT__
299f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
300f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
301f3fbd535SBarry Smith {
302f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
303f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
304f3fbd535SBarry Smith   PetscErrorCode ierr;
305f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
306f3fbd535SBarry Smith 
307f3fbd535SBarry Smith   PetscFunctionBegin;
308e1d8e5deSBarry Smith 
309e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
310298cc208SBarry Smith   for (i=0; i<levels; i++) {
311e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
312e1d8e5deSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
313298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
314e1d8e5deSBarry Smith     }
315e1d8e5deSBarry Smith   }
316e1d8e5deSBarry Smith 
317f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
318f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
31931567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
320f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
32131567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
322f3fbd535SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
323f3fbd535SBarry Smith     }
324f3fbd535SBarry Smith   }
32531567311SBarry Smith   else if (mg->am == PC_MG_ADDITIVE) {
32631567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
327f3fbd535SBarry Smith   }
32831567311SBarry Smith   else if (mg->am == PC_MG_KASKADE) {
32931567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
330f3fbd535SBarry Smith   }
331f3fbd535SBarry Smith   else {
332f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
333f3fbd535SBarry Smith   }
334f3fbd535SBarry Smith   PetscFunctionReturn(0);
335f3fbd535SBarry Smith }
336f3fbd535SBarry Smith 
337f3fbd535SBarry Smith 
338f3fbd535SBarry Smith #undef __FUNCT__
339f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
340f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc)
341f3fbd535SBarry Smith {
342f3fbd535SBarry Smith   PetscErrorCode ierr;
343f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
344ace3abfcSBarry Smith   PetscBool      flg;
345f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
346f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
347f3fbd535SBarry Smith   PCMGType       mgtype;
348f3fbd535SBarry Smith   PCMGCycleType  mgctype;
349f3fbd535SBarry Smith 
350f3fbd535SBarry Smith   PetscFunctionBegin;
351f3fbd535SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
3523aeaf226SBarry Smith     if (!mg->levels) {
353f3fbd535SBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
354eb3f98d2SBarry Smith       if (!flg && pc->dm) {
355eb3f98d2SBarry Smith         ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
356eb3f98d2SBarry Smith         levels++;
3573aeaf226SBarry Smith         mg->usedmfornumberoflevels = PETSC_TRUE;
358eb3f98d2SBarry Smith       }
359f3fbd535SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
360f3fbd535SBarry Smith     }
3613aeaf226SBarry Smith     mglevels = mg->levels;
3623aeaf226SBarry Smith 
363f3fbd535SBarry Smith     mgctype = (PCMGCycleType) mglevels[0]->cycles;
364f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
365f3fbd535SBarry Smith     if (flg) {
366f3fbd535SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
367f3fbd535SBarry Smith     };
368f3fbd535SBarry Smith     flg  = PETSC_FALSE;
369acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
370f3fbd535SBarry Smith     if (flg) {
371f3fbd535SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
372f3fbd535SBarry Smith     }
373f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
374f3fbd535SBarry Smith     if (flg) {
375f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
376f3fbd535SBarry Smith     }
377f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
378f3fbd535SBarry Smith     if (flg) {
379f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
380f3fbd535SBarry Smith     }
38131567311SBarry Smith     mgtype = mg->am;
382f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
383f3fbd535SBarry Smith     if (flg) {
384f3fbd535SBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
385f3fbd535SBarry Smith     }
38631567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
38731567311SBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
388f3fbd535SBarry Smith       if (flg) {
389f3fbd535SBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
390f3fbd535SBarry Smith       }
391f3fbd535SBarry Smith     }
392f3fbd535SBarry Smith     flg  = PETSC_FALSE;
393acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
394f3fbd535SBarry Smith     if (flg) {
395f3fbd535SBarry Smith       PetscInt i;
396f3fbd535SBarry Smith       char     eventname[128];
39763e6d426SJed Brown       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
398f3fbd535SBarry Smith       levels = mglevels[0]->levels;
399f3fbd535SBarry Smith       for (i=0; i<levels; i++) {
400f3fbd535SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
40163e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
402f3fbd535SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
40363e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
404f3fbd535SBarry Smith         if (i) {
405f3fbd535SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
40663e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
407f3fbd535SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
40863e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
409f3fbd535SBarry Smith         }
410f3fbd535SBarry Smith       }
411f3fbd535SBarry Smith     }
412f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
413f3fbd535SBarry Smith   PetscFunctionReturn(0);
414f3fbd535SBarry Smith }
415f3fbd535SBarry Smith 
416f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
417f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
418f3fbd535SBarry Smith 
419f3fbd535SBarry Smith #undef __FUNCT__
420f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
421f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
422f3fbd535SBarry Smith {
423f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
424f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
425f3fbd535SBarry Smith   PetscErrorCode ierr;
426f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
427ace3abfcSBarry Smith   PetscBool      iascii;
428f3fbd535SBarry Smith 
429f3fbd535SBarry Smith   PetscFunctionBegin;
4302692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
431f3fbd535SBarry Smith   if (iascii) {
43231567311SBarry 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);
43331567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
43431567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
435f3fbd535SBarry Smith     }
436218a07d4SBarry Smith     if (mg->galerkin) {
437f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4384f66f45eSBarry Smith     } else {
4394f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
440f3fbd535SBarry Smith     }
441f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
442f3fbd535SBarry Smith       if (!i) {
4432d19a89fSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
444f3fbd535SBarry Smith       } else {
44531567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
446f3fbd535SBarry Smith       }
447f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
448f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
449f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
450f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
451f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
452f3fbd535SBarry Smith       } else if (i){
45331567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr);
454f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
455f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
456f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
457f3fbd535SBarry Smith       }
458f3fbd535SBarry Smith     }
459f3fbd535SBarry Smith   } else {
46065e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
461f3fbd535SBarry Smith   }
462f3fbd535SBarry Smith   PetscFunctionReturn(0);
463f3fbd535SBarry Smith }
464f3fbd535SBarry Smith 
465cab2e9ccSBarry Smith #include <private/dmimpl.h>
466cab2e9ccSBarry Smith #include <private/kspimpl.h>
467cab2e9ccSBarry Smith 
468f3fbd535SBarry Smith /*
469f3fbd535SBarry Smith     Calls setup for the KSP on each level
470f3fbd535SBarry Smith */
471f3fbd535SBarry Smith #undef __FUNCT__
472f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
473f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
474f3fbd535SBarry Smith {
475f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
476f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
477f3fbd535SBarry Smith   PetscErrorCode          ierr;
478f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
479f3fbd535SBarry Smith   PC                      cpc,mpc;
480fd1303e9SJungho Lee   PetscBool               preonly,lu,redundant,cholesky,svd,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset;
481f3fbd535SBarry Smith   PetscViewerASCIIMonitor ascii;
482f3fbd535SBarry Smith   PetscViewer             viewer = PETSC_NULL;
483f3fbd535SBarry Smith   MPI_Comm                comm;
484f3fbd535SBarry Smith   Mat                     dA,dB;
485f3fbd535SBarry Smith   MatStructure            uflag;
486f3fbd535SBarry Smith   Vec                     tvec;
487218a07d4SBarry Smith   DM                      *dms;
488f3fbd535SBarry Smith 
489f3fbd535SBarry Smith   PetscFunctionBegin;
4903aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
4913aeaf226SBarry Smith     PetscInt levels;
4923aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
4933aeaf226SBarry Smith     levels++;
4943aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
4953aeaf226SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
4963aeaf226SBarry Smith       n    = levels;
4973aeaf226SBarry Smith       ierr = PCSetFromOptions(pc);CHKERRQ(ierr);  /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
4983aeaf226SBarry Smith       mglevels =  mg->levels;
4993aeaf226SBarry Smith     }
5003aeaf226SBarry Smith   }
5013aeaf226SBarry Smith 
502f3fbd535SBarry Smith 
503f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
504f3fbd535SBarry Smith   /* so use those from global PC */
505f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
506f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
507f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
508f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr);
5091338a6b9SJed Brown   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) {
510f3fbd535SBarry 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);
511f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
512f3fbd535SBarry Smith   }
513f3fbd535SBarry Smith 
5142d2b81a6SBarry Smith   if (pc->dm && !pc->setupcalled) {
5152d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
5162e499ae9SBarry Smith     Mat p;
51773250ac0SBarry Smith     Vec rscale;
518218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
519218a07d4SBarry Smith     dms[n-1] = pc->dm;
520218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
521218a07d4SBarry Smith       ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr);
5223c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
5233c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
52411629dbeSBarry Smith       ierr = DMSetFunction(dms[i],0);
52511629dbeSBarry Smith       ierr = DMSetInitialGuess(dms[i],0);
52624c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
52773250ac0SBarry Smith 	ierr = DMGetInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
5282d2b81a6SBarry Smith 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
529*00ac8be1SBarry Smith 	if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
53073250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
5316bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
532218a07d4SBarry Smith       }
53324c3aa18SBarry Smith     }
5342d2b81a6SBarry Smith 
535218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
5366bf464f9SBarry Smith       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
537218a07d4SBarry Smith     }
5382d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
539cab2e9ccSBarry Smith 
540cab2e9ccSBarry Smith     /* finest smoother also gets DM but it is not active */
541cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
542cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
543218a07d4SBarry Smith   }
544218a07d4SBarry Smith 
545218a07d4SBarry Smith   if (mg->galerkin) {
546f3fbd535SBarry Smith     Mat B;
547218a07d4SBarry Smith     mg->galerkinused = PETSC_TRUE;
548f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
549f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
550f3fbd535SBarry Smith     if (!pc->setupcalled) {
551f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
552f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
553f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
554f3fbd535SBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
555f3fbd535SBarry Smith         dB   = B;
556f3fbd535SBarry Smith       }
557cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
558f3fbd535SBarry Smith     } else {
559f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
560f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
561f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
562f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
563f3fbd535SBarry Smith         dB   = B;
564f3fbd535SBarry Smith       }
565f3fbd535SBarry Smith     }
566cab2e9ccSBarry Smith   } else if (pc->dm && pc->dm->x) {
567cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
568cab2e9ccSBarry Smith     for (i=n-2;i>-1; i--) {
569cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
570cab2e9ccSBarry Smith         Vec *vecs;
571cab2e9ccSBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr);
572cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
573cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
574cab2e9ccSBarry Smith       }
575cab2e9ccSBarry Smith       ierr = MatRestrict(mglevels[i+1]->interpolate,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
57673250ac0SBarry Smith       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,mglevels[i+1]->rscale);CHKERRQ(ierr);
577cab2e9ccSBarry Smith     }
578f3fbd535SBarry Smith   }
579f3fbd535SBarry Smith 
580f3fbd535SBarry Smith   if (!pc->setupcalled) {
581acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr);
582f3fbd535SBarry Smith 
583f3fbd535SBarry Smith     for (i=0; i<n; i++) {
584f3fbd535SBarry Smith       if (monitor) {
585f3fbd535SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr);
586f3fbd535SBarry Smith         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
587c2efdce3SBarry Smith         ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
588f3fbd535SBarry Smith       }
589f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
590f3fbd535SBarry Smith     }
591f3fbd535SBarry Smith     for (i=1; i<n; i++) {
592f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
593f3fbd535SBarry Smith         if (monitor) {
594f3fbd535SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr);
595f3fbd535SBarry Smith           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
596c2efdce3SBarry Smith           ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
597f3fbd535SBarry Smith         }
598f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
599f3fbd535SBarry Smith       }
600f3fbd535SBarry Smith     }
601f3fbd535SBarry Smith     for (i=1; i<n; i++) {
602f3fbd535SBarry Smith       if (mglevels[i]->restrct && !mglevels[i]->interpolate) {
603f3fbd535SBarry Smith         ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr);
604f3fbd535SBarry Smith       }
605f3fbd535SBarry Smith       if (!mglevels[i]->restrct && mglevels[i]->interpolate) {
606f3fbd535SBarry Smith         ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr);
607f3fbd535SBarry Smith       }
608f3fbd535SBarry Smith #if defined(PETSC_USE_DEBUG)
609f3fbd535SBarry Smith       if (!mglevels[i]->restrct || !mglevels[i]->interpolate) {
61065e19b50SBarry Smith         SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
611f3fbd535SBarry Smith       }
612f3fbd535SBarry Smith #endif
613f3fbd535SBarry Smith     }
614f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
615f3fbd535SBarry Smith       if (!mglevels[i]->b) {
616f3fbd535SBarry Smith         Vec *vec;
617f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
618f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
6196bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
620f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
621f3fbd535SBarry Smith       }
622f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
623f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
624f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
6256bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
626f3fbd535SBarry Smith       }
627f3fbd535SBarry Smith       if (!mglevels[i]->x) {
628f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
629f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
6306bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
631f3fbd535SBarry Smith       }
632f3fbd535SBarry Smith     }
633f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
634f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
635f3fbd535SBarry Smith       Vec *vec;
636f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
637f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
6386bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
639f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
640f3fbd535SBarry Smith     }
641f3fbd535SBarry Smith   }
642f3fbd535SBarry Smith 
643f3fbd535SBarry Smith 
644f3fbd535SBarry Smith   for (i=1; i<n; i++) {
645f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
646f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
647f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
648f3fbd535SBarry Smith     }
64963e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
650f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
65163e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
652d42688cbSBarry Smith     if (!mglevels[i]->residual) {
653d42688cbSBarry Smith       Mat mat;
654d42688cbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
655d42688cbSBarry Smith       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
656d42688cbSBarry Smith     }
657f3fbd535SBarry Smith   }
658f3fbd535SBarry Smith   for (i=1; i<n; i++) {
659f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
660f3fbd535SBarry Smith       Mat          downmat,downpmat;
661f3fbd535SBarry Smith       MatStructure matflag;
662ace3abfcSBarry Smith       PetscBool    opsset;
663f3fbd535SBarry Smith 
664f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
665f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
666f3fbd535SBarry Smith       if (!opsset) {
667f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
668f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
669f3fbd535SBarry Smith       }
670f3fbd535SBarry Smith 
671f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
67263e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
673f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
67463e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
675f3fbd535SBarry Smith     }
676f3fbd535SBarry Smith   }
677f3fbd535SBarry Smith 
678f3fbd535SBarry Smith   /*
679f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
680f3fbd535SBarry Smith   */
681f3fbd535SBarry Smith   ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
682f3fbd535SBarry Smith   if (preonly) {
683f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
684f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
685f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
686fd1303e9SJungho Lee     ierr = PetscTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
687fd1303e9SJungho Lee     if (!lu && !redundant && !cholesky && !svd) {
688f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
689f3fbd535SBarry Smith     }
690f3fbd535SBarry Smith   }
691f3fbd535SBarry Smith 
692f3fbd535SBarry Smith   if (!pc->setupcalled) {
693f3fbd535SBarry Smith     if (monitor) {
694f3fbd535SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr);
695f3fbd535SBarry Smith       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
696c2efdce3SBarry Smith       ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
697f3fbd535SBarry Smith     }
698f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
699f3fbd535SBarry Smith   }
700f3fbd535SBarry Smith 
70163e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
702f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
70363e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
704f3fbd535SBarry Smith 
705f3fbd535SBarry Smith   /*
706f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
707e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
708f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
709f3fbd535SBarry Smith 
710f3fbd535SBarry Smith    Only support one or the other at the same time.
711f3fbd535SBarry Smith   */
712f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
713acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
714f3fbd535SBarry Smith   if (dump) {
715f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
716f3fbd535SBarry Smith   }
717f3fbd535SBarry Smith   dump = PETSC_FALSE;
718f3fbd535SBarry Smith #endif
719acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
720f3fbd535SBarry Smith   if (dump) {
721f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
722f3fbd535SBarry Smith   }
723f3fbd535SBarry Smith 
724f3fbd535SBarry Smith   if (viewer) {
725f3fbd535SBarry Smith     for (i=1; i<n; i++) {
726f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
727f3fbd535SBarry Smith     }
728f3fbd535SBarry Smith     for (i=0; i<n; i++) {
729f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
730f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
731f3fbd535SBarry Smith     }
732f3fbd535SBarry Smith   }
733f3fbd535SBarry Smith   PetscFunctionReturn(0);
734f3fbd535SBarry Smith }
735f3fbd535SBarry Smith 
736f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
737f3fbd535SBarry Smith 
738f3fbd535SBarry Smith #undef __FUNCT__
7399dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
7404b9ad928SBarry Smith /*@
74197177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
7424b9ad928SBarry Smith 
7434b9ad928SBarry Smith    Not Collective
7444b9ad928SBarry Smith 
7454b9ad928SBarry Smith    Input Parameter:
7464b9ad928SBarry Smith .  pc - the preconditioner context
7474b9ad928SBarry Smith 
7484b9ad928SBarry Smith    Output parameter:
7494b9ad928SBarry Smith .  levels - the number of levels
7504b9ad928SBarry Smith 
7514b9ad928SBarry Smith    Level: advanced
7524b9ad928SBarry Smith 
7534b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
7544b9ad928SBarry Smith 
75597177400SBarry Smith .seealso: PCMGSetLevels()
7564b9ad928SBarry Smith @*/
7577087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
7584b9ad928SBarry Smith {
759f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
7604b9ad928SBarry Smith 
7614b9ad928SBarry Smith   PetscFunctionBegin;
7620700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
7634482741eSBarry Smith   PetscValidIntPointer(levels,2);
764f3fbd535SBarry Smith   *levels = mg->nlevels;
7654b9ad928SBarry Smith   PetscFunctionReturn(0);
7664b9ad928SBarry Smith }
7674b9ad928SBarry Smith 
7684b9ad928SBarry Smith #undef __FUNCT__
7699dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
7704b9ad928SBarry Smith /*@
77197177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
7724b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
7734b9ad928SBarry Smith 
774ad4df100SBarry Smith    Logically Collective on PC
7754b9ad928SBarry Smith 
7764b9ad928SBarry Smith    Input Parameters:
7774b9ad928SBarry Smith +  pc - the preconditioner context
7789dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
7799dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
7804b9ad928SBarry Smith 
7814b9ad928SBarry Smith    Options Database Key:
7824b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
7834b9ad928SBarry Smith    additive, full, kaskade
7844b9ad928SBarry Smith 
7854b9ad928SBarry Smith    Level: advanced
7864b9ad928SBarry Smith 
7874b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
7884b9ad928SBarry Smith 
78997177400SBarry Smith .seealso: PCMGSetLevels()
7904b9ad928SBarry Smith @*/
7917087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
7924b9ad928SBarry Smith {
793f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
7944b9ad928SBarry Smith 
7954b9ad928SBarry Smith   PetscFunctionBegin;
7960700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
797c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
79831567311SBarry Smith   mg->am = form;
7999dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
8004b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
8014b9ad928SBarry Smith   PetscFunctionReturn(0);
8024b9ad928SBarry Smith }
8034b9ad928SBarry Smith 
8044b9ad928SBarry Smith #undef __FUNCT__
8050d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
8064b9ad928SBarry Smith /*@
8070d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
8084b9ad928SBarry Smith    complicated cycling.
8094b9ad928SBarry Smith 
810ad4df100SBarry Smith    Logically Collective on PC
8114b9ad928SBarry Smith 
8124b9ad928SBarry Smith    Input Parameters:
813c2be2410SBarry Smith +  pc - the multigrid context
8140d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
8154b9ad928SBarry Smith 
8164b9ad928SBarry Smith    Options Database Key:
8170d353602SBarry Smith $  -pc_mg_cycle_type v or w
8184b9ad928SBarry Smith 
8194b9ad928SBarry Smith    Level: advanced
8204b9ad928SBarry Smith 
8214b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8224b9ad928SBarry Smith 
8230d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
8244b9ad928SBarry Smith @*/
8257087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
8264b9ad928SBarry Smith {
827f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
828f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
82979416396SBarry Smith   PetscInt     i,levels;
8304b9ad928SBarry Smith 
8314b9ad928SBarry Smith   PetscFunctionBegin;
8320700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
833e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
834c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
835f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8364b9ad928SBarry Smith 
8374b9ad928SBarry Smith   for (i=0; i<levels; i++) {
838f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
8394b9ad928SBarry Smith   }
8404b9ad928SBarry Smith   PetscFunctionReturn(0);
8414b9ad928SBarry Smith }
8424b9ad928SBarry Smith 
8434b9ad928SBarry Smith #undef __FUNCT__
8448cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
8458cc2d5dfSBarry Smith /*@
8468cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
8478cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
8488cc2d5dfSBarry Smith 
849ad4df100SBarry Smith    Logically Collective on PC
8508cc2d5dfSBarry Smith 
8518cc2d5dfSBarry Smith    Input Parameters:
8528cc2d5dfSBarry Smith +  pc - the multigrid context
8538cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
8548cc2d5dfSBarry Smith 
8558cc2d5dfSBarry Smith    Options Database Key:
8568cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
8578cc2d5dfSBarry Smith 
8588cc2d5dfSBarry Smith    Level: advanced
8598cc2d5dfSBarry Smith 
8608cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
8618cc2d5dfSBarry Smith 
8628cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8638cc2d5dfSBarry Smith 
8648cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
8658cc2d5dfSBarry Smith @*/
8667087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
8678cc2d5dfSBarry Smith {
868f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
869f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
8708cc2d5dfSBarry Smith   PetscInt     i,levels;
8718cc2d5dfSBarry Smith 
8728cc2d5dfSBarry Smith   PetscFunctionBegin;
8730700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
874e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
875c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
876f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8778cc2d5dfSBarry Smith 
8788cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
87931567311SBarry Smith     mg->cyclesperpcapply  = n;
8808cc2d5dfSBarry Smith   }
8818cc2d5dfSBarry Smith   PetscFunctionReturn(0);
8828cc2d5dfSBarry Smith }
8838cc2d5dfSBarry Smith 
8848cc2d5dfSBarry Smith #undef __FUNCT__
8859dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
886c2be2410SBarry Smith /*@
88797177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
888c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
889c2be2410SBarry Smith 
890ad4df100SBarry Smith    Logically Collective on PC
891c2be2410SBarry Smith 
892c2be2410SBarry Smith    Input Parameters:
8933fc8bf9cSBarry Smith .  pc - the multigrid context
894c2be2410SBarry Smith 
895c2be2410SBarry Smith    Options Database Key:
896c2be2410SBarry Smith $  -pc_mg_galerkin
897c2be2410SBarry Smith 
898c2be2410SBarry Smith    Level: intermediate
899c2be2410SBarry Smith 
900c2be2410SBarry Smith .keywords: MG, set, Galerkin
901c2be2410SBarry Smith 
9023fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
9033fc8bf9cSBarry Smith 
904c2be2410SBarry Smith @*/
9057087cfbeSBarry Smith PetscErrorCode  PCMGSetGalerkin(PC pc)
906c2be2410SBarry Smith {
907f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
908c2be2410SBarry Smith 
909c2be2410SBarry Smith   PetscFunctionBegin;
9100700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
911218a07d4SBarry Smith   mg->galerkin = PETSC_TRUE;
912c2be2410SBarry Smith   PetscFunctionReturn(0);
913c2be2410SBarry Smith }
914c2be2410SBarry Smith 
915c2be2410SBarry Smith #undef __FUNCT__
9163fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
9173fc8bf9cSBarry Smith /*@
9183fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
9193fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
9203fc8bf9cSBarry Smith 
9213fc8bf9cSBarry Smith    Not Collective
9223fc8bf9cSBarry Smith 
9233fc8bf9cSBarry Smith    Input Parameter:
9243fc8bf9cSBarry Smith .  pc - the multigrid context
9253fc8bf9cSBarry Smith 
9263fc8bf9cSBarry Smith    Output Parameter:
9273fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
9283fc8bf9cSBarry Smith 
9293fc8bf9cSBarry Smith    Options Database Key:
9303fc8bf9cSBarry Smith $  -pc_mg_galerkin
9313fc8bf9cSBarry Smith 
9323fc8bf9cSBarry Smith    Level: intermediate
9333fc8bf9cSBarry Smith 
9343fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
9353fc8bf9cSBarry Smith 
9363fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
9373fc8bf9cSBarry Smith 
9383fc8bf9cSBarry Smith @*/
9397087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
9403fc8bf9cSBarry Smith {
941f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
9423fc8bf9cSBarry Smith 
9433fc8bf9cSBarry Smith   PetscFunctionBegin;
9440700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
945218a07d4SBarry Smith   *galerkin = mg->galerkin;
9463fc8bf9cSBarry Smith   PetscFunctionReturn(0);
9473fc8bf9cSBarry Smith }
9483fc8bf9cSBarry Smith 
9493fc8bf9cSBarry Smith #undef __FUNCT__
9509dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
9514b9ad928SBarry Smith /*@
95297177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
95397177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
9544b9ad928SBarry Smith    pre-smoothing steps on different levels.
9554b9ad928SBarry Smith 
956ad4df100SBarry Smith    Logically Collective on PC
9574b9ad928SBarry Smith 
9584b9ad928SBarry Smith    Input Parameters:
9594b9ad928SBarry Smith +  mg - the multigrid context
9604b9ad928SBarry Smith -  n - the number of smoothing steps
9614b9ad928SBarry Smith 
9624b9ad928SBarry Smith    Options Database Key:
9634b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
9644b9ad928SBarry Smith 
9654b9ad928SBarry Smith    Level: advanced
9664b9ad928SBarry Smith 
9674b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
9684b9ad928SBarry Smith 
96997177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
9704b9ad928SBarry Smith @*/
9717087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
9724b9ad928SBarry Smith {
973f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
974f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9756849ba73SBarry Smith   PetscErrorCode ierr;
97679416396SBarry Smith   PetscInt       i,levels;
9774b9ad928SBarry Smith 
9784b9ad928SBarry Smith   PetscFunctionBegin;
9790700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
980e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
981c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
982f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9834b9ad928SBarry Smith 
984b05257ddSBarry Smith   for (i=1; i<levels; i++) {
9854b9ad928SBarry Smith     /* make sure smoother up and down are different */
98697177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
987f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
98831567311SBarry Smith     mg->default_smoothd = n;
9894b9ad928SBarry Smith   }
9904b9ad928SBarry Smith   PetscFunctionReturn(0);
9914b9ad928SBarry Smith }
9924b9ad928SBarry Smith 
9934b9ad928SBarry Smith #undef __FUNCT__
9949dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
9954b9ad928SBarry Smith /*@
99697177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
99797177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
9984b9ad928SBarry Smith    post-smoothing steps on different levels.
9994b9ad928SBarry Smith 
1000ad4df100SBarry Smith    Logically Collective on PC
10014b9ad928SBarry Smith 
10024b9ad928SBarry Smith    Input Parameters:
10034b9ad928SBarry Smith +  mg - the multigrid context
10044b9ad928SBarry Smith -  n - the number of smoothing steps
10054b9ad928SBarry Smith 
10064b9ad928SBarry Smith    Options Database Key:
10074b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
10084b9ad928SBarry Smith 
10094b9ad928SBarry Smith    Level: advanced
10104b9ad928SBarry Smith 
10114b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
1012a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
10134b9ad928SBarry Smith 
10144b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
10154b9ad928SBarry Smith 
101697177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
10174b9ad928SBarry Smith @*/
10187087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
10194b9ad928SBarry Smith {
1020f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1021f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10226849ba73SBarry Smith   PetscErrorCode ierr;
102379416396SBarry Smith   PetscInt       i,levels;
10244b9ad928SBarry Smith 
10254b9ad928SBarry Smith   PetscFunctionBegin;
10260700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1027e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1028c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1029f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10304b9ad928SBarry Smith 
10314b9ad928SBarry Smith   for (i=1; i<levels; i++) {
10324b9ad928SBarry Smith     /* make sure smoother up and down are different */
103397177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1034f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
103531567311SBarry Smith     mg->default_smoothu = n;
10364b9ad928SBarry Smith   }
10374b9ad928SBarry Smith   PetscFunctionReturn(0);
10384b9ad928SBarry Smith }
10394b9ad928SBarry Smith 
10404b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
10414b9ad928SBarry Smith 
10423b09bd56SBarry Smith /*MC
1043ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
10443b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
10453b09bd56SBarry Smith 
10463b09bd56SBarry Smith    Options Database Keys:
10473b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
10480d353602SBarry Smith .  -pc_mg_cycles v or w
104979416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
10503b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
10513b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
10523b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
10533b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
105468eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
10553b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1056e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
10573b09bd56SBarry Smith 
105824c3aa18SBarry 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
10593b09bd56SBarry Smith 
10603b09bd56SBarry Smith    Level: intermediate
10613b09bd56SBarry Smith 
10628f87f92bSBarry Smith    Concepts: multigrid/multilevel
10633b09bd56SBarry Smith 
106424c3aa18SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC
10650d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
106697177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
106797177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
10680d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
10693b09bd56SBarry Smith M*/
10703b09bd56SBarry Smith 
10714b9ad928SBarry Smith EXTERN_C_BEGIN
10724b9ad928SBarry Smith #undef __FUNCT__
10734b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
10747087cfbeSBarry Smith PetscErrorCode  PCCreate_MG(PC pc)
10754b9ad928SBarry Smith {
1076f3fbd535SBarry Smith   PC_MG          *mg;
1077f3fbd535SBarry Smith   PetscErrorCode ierr;
1078f3fbd535SBarry Smith 
10794b9ad928SBarry Smith   PetscFunctionBegin;
1080f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1081f3fbd535SBarry Smith   pc->data    = (void*)mg;
1082f3fbd535SBarry Smith   mg->nlevels = -1;
1083f3fbd535SBarry Smith 
10844b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
10854b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1086a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
10874b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
10884b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
10894b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
10904b9ad928SBarry Smith   PetscFunctionReturn(0);
10914b9ad928SBarry Smith }
10924b9ad928SBarry Smith EXTERN_C_END
1093