xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision c88c5224601c4a9900c397283ffbc1685ca8c72c)
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;
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 
205f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
206f3fbd535SBarry Smith     ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr);
207f3fbd535SBarry Smith     mglevels[i]->level           = i;
208f3fbd535SBarry Smith     mglevels[i]->levels          = levels;
209f3fbd535SBarry Smith     mglevels[i]->cycles          = PC_MG_CYCLE_V;
21031567311SBarry Smith     mg->default_smoothu = 1;
21131567311SBarry Smith     mg->default_smoothd = 1;
21263e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
21363e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
21463e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
21563e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
216f3fbd535SBarry Smith 
217f3fbd535SBarry Smith     if (comms) comm = comms[i];
218f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
219f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
22031567311SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
221f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
222f3fbd535SBarry Smith 
223f3fbd535SBarry Smith     /* do special stuff for coarse grid */
224f3fbd535SBarry Smith     if (!i && levels > 1) {
225f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
226f3fbd535SBarry Smith 
227f3fbd535SBarry Smith       /* coarse solve is (redundant) LU by default */
228f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
229f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
230f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
231f3fbd535SBarry Smith       if (size > 1) {
232f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
233f3fbd535SBarry Smith       } else {
234f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
235f3fbd535SBarry Smith       }
236f3fbd535SBarry Smith 
237f3fbd535SBarry Smith     } else {
238f3fbd535SBarry Smith       char tprefix[128];
239f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
240f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
241f3fbd535SBarry Smith     }
242f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr);
243f3fbd535SBarry Smith     mglevels[i]->smoothu    = mglevels[i]->smoothd;
24431567311SBarry Smith     mg->rtol                = 0.0;
24531567311SBarry Smith     mg->abstol              = 0.0;
24631567311SBarry Smith     mg->dtol                = 0.0;
24731567311SBarry Smith     mg->ttol                = 0.0;
24831567311SBarry Smith     mg->cyclesperpcapply    = 1;
249f3fbd535SBarry Smith   }
25031567311SBarry Smith   mg->am          = PC_MG_MULTIPLICATIVE;
251f3fbd535SBarry Smith   mg->levels      = mglevels;
2524b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
2534b9ad928SBarry Smith   PetscFunctionReturn(0);
2544b9ad928SBarry Smith }
2554b9ad928SBarry Smith 
256c07bf074SBarry Smith 
257c07bf074SBarry Smith #undef __FUNCT__
258c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
259c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
260c07bf074SBarry Smith {
261c07bf074SBarry Smith   PetscErrorCode ierr;
262a06653b4SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
263a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
264a06653b4SBarry Smith   PetscInt       i,n;
265c07bf074SBarry Smith 
266c07bf074SBarry Smith   PetscFunctionBegin;
267a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
268a06653b4SBarry Smith   if (mglevels) {
269a06653b4SBarry Smith     n = mglevels[0]->levels;
270a06653b4SBarry Smith     for (i=0; i<n; i++) {
271a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2726bf464f9SBarry Smith 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
273a06653b4SBarry Smith       }
2746bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
275a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
276a06653b4SBarry Smith     }
277a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
278a06653b4SBarry Smith   }
279c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
280f3fbd535SBarry Smith   PetscFunctionReturn(0);
281f3fbd535SBarry Smith }
282f3fbd535SBarry Smith 
283f3fbd535SBarry Smith 
284f3fbd535SBarry Smith 
28509573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
28609573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
28709573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
288f3fbd535SBarry Smith 
289f3fbd535SBarry Smith /*
290f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
291f3fbd535SBarry Smith              or full cycle of multigrid.
292f3fbd535SBarry Smith 
293f3fbd535SBarry Smith   Note:
294f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
295f3fbd535SBarry Smith */
296f3fbd535SBarry Smith #undef __FUNCT__
297f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
298f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
299f3fbd535SBarry Smith {
300f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
301f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
302f3fbd535SBarry Smith   PetscErrorCode ierr;
303f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
304f3fbd535SBarry Smith 
305f3fbd535SBarry Smith   PetscFunctionBegin;
306e1d8e5deSBarry Smith 
307e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
308298cc208SBarry Smith   for (i=0; i<levels; i++) {
309e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
310e1d8e5deSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
311298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
312e1d8e5deSBarry Smith     }
313e1d8e5deSBarry Smith   }
314e1d8e5deSBarry Smith 
315f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
316f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
31731567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
318f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
31931567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
320f3fbd535SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
321f3fbd535SBarry Smith     }
322f3fbd535SBarry Smith   }
32331567311SBarry Smith   else if (mg->am == PC_MG_ADDITIVE) {
32431567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
325f3fbd535SBarry Smith   }
32631567311SBarry Smith   else if (mg->am == PC_MG_KASKADE) {
32731567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
328f3fbd535SBarry Smith   }
329f3fbd535SBarry Smith   else {
330f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
331f3fbd535SBarry Smith   }
332f3fbd535SBarry Smith   PetscFunctionReturn(0);
333f3fbd535SBarry Smith }
334f3fbd535SBarry Smith 
335f3fbd535SBarry Smith 
336f3fbd535SBarry Smith #undef __FUNCT__
337f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
338f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc)
339f3fbd535SBarry Smith {
340f3fbd535SBarry Smith   PetscErrorCode ierr;
341f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
342c91913d3SJed Brown   PetscBool      flg,set;
343f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
344f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
345f3fbd535SBarry Smith   PCMGType       mgtype;
346f3fbd535SBarry Smith   PCMGCycleType  mgctype;
347f3fbd535SBarry Smith 
348f3fbd535SBarry Smith   PetscFunctionBegin;
349f3fbd535SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
3503aeaf226SBarry Smith     if (!mg->levels) {
351f3fbd535SBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
352eb3f98d2SBarry Smith       if (!flg && pc->dm) {
353eb3f98d2SBarry Smith         ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
354eb3f98d2SBarry Smith         levels++;
3553aeaf226SBarry Smith         mg->usedmfornumberoflevels = PETSC_TRUE;
356eb3f98d2SBarry Smith       }
357f3fbd535SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
358f3fbd535SBarry Smith     }
3593aeaf226SBarry Smith     mglevels = mg->levels;
3603aeaf226SBarry Smith 
361f3fbd535SBarry Smith     mgctype = (PCMGCycleType) mglevels[0]->cycles;
362f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
363f3fbd535SBarry Smith     if (flg) {
364f3fbd535SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
365f3fbd535SBarry Smith     };
366f3fbd535SBarry Smith     flg  = PETSC_FALSE;
367c91913d3SJed Brown     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr);
368c91913d3SJed Brown     if (set) {
369c91913d3SJed Brown       ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr);
370f3fbd535SBarry Smith     }
371f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
372f3fbd535SBarry Smith     if (flg) {
373f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
374f3fbd535SBarry Smith     }
375f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
376f3fbd535SBarry Smith     if (flg) {
377f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
378f3fbd535SBarry Smith     }
37931567311SBarry Smith     mgtype = mg->am;
380f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
381f3fbd535SBarry Smith     if (flg) {
382f3fbd535SBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
383f3fbd535SBarry Smith     }
38431567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
38531567311SBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
386f3fbd535SBarry Smith       if (flg) {
387f3fbd535SBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
388f3fbd535SBarry Smith       }
389f3fbd535SBarry Smith     }
390f3fbd535SBarry Smith     flg  = PETSC_FALSE;
391acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
392f3fbd535SBarry Smith     if (flg) {
393f3fbd535SBarry Smith       PetscInt i;
394f3fbd535SBarry Smith       char     eventname[128];
39563e6d426SJed Brown       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
396f3fbd535SBarry Smith       levels = mglevels[0]->levels;
397f3fbd535SBarry Smith       for (i=0; i<levels; i++) {
398f3fbd535SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
39963e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
400f3fbd535SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
40163e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
402f3fbd535SBarry Smith         if (i) {
403f3fbd535SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
40463e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
405f3fbd535SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
40663e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
407f3fbd535SBarry Smith         }
408f3fbd535SBarry Smith       }
409f3fbd535SBarry Smith     }
410f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
411f3fbd535SBarry Smith   PetscFunctionReturn(0);
412f3fbd535SBarry Smith }
413f3fbd535SBarry Smith 
414f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
415f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
416f3fbd535SBarry Smith 
417f3fbd535SBarry Smith #undef __FUNCT__
418f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
419f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
420f3fbd535SBarry Smith {
421f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
422f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
423f3fbd535SBarry Smith   PetscErrorCode ierr;
424f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
425ace3abfcSBarry Smith   PetscBool      iascii;
426f3fbd535SBarry Smith 
427f3fbd535SBarry Smith   PetscFunctionBegin;
4282692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
429f3fbd535SBarry Smith   if (iascii) {
43031567311SBarry 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);
43131567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
43231567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
433f3fbd535SBarry Smith     }
434218a07d4SBarry Smith     if (mg->galerkin) {
435f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4364f66f45eSBarry Smith     } else {
4374f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
438f3fbd535SBarry Smith     }
439f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
440f3fbd535SBarry Smith       if (!i) {
44189cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
442f3fbd535SBarry Smith       } else {
44389cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
444f3fbd535SBarry Smith       }
445f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
446f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
447f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
448f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
449f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
450f3fbd535SBarry Smith       } else if (i){
45189cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
452f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
453f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
454f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
455f3fbd535SBarry Smith       }
456f3fbd535SBarry Smith     }
457649052a6SBarry Smith   } else SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
458f3fbd535SBarry Smith   PetscFunctionReturn(0);
459f3fbd535SBarry Smith }
460f3fbd535SBarry Smith 
461cab2e9ccSBarry Smith #include <private/dmimpl.h>
462cab2e9ccSBarry Smith #include <private/kspimpl.h>
463cab2e9ccSBarry Smith 
464f3fbd535SBarry Smith /*
465f3fbd535SBarry Smith     Calls setup for the KSP on each level
466f3fbd535SBarry Smith */
467f3fbd535SBarry Smith #undef __FUNCT__
468f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
469f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
470f3fbd535SBarry Smith {
471f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
472f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
473f3fbd535SBarry Smith   PetscErrorCode          ierr;
474f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
475f3fbd535SBarry Smith   PC                      cpc,mpc;
476649052a6SBarry Smith   PetscBool               preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset;
477f3fbd535SBarry Smith   Mat                     dA,dB;
478f3fbd535SBarry Smith   MatStructure            uflag;
479f3fbd535SBarry Smith   Vec                     tvec;
480218a07d4SBarry Smith   DM                      *dms;
481649052a6SBarry Smith   PetscViewer             viewer = 0;
482f3fbd535SBarry Smith 
483f3fbd535SBarry Smith   PetscFunctionBegin;
4843aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
4853aeaf226SBarry Smith     PetscInt levels;
4863aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
4873aeaf226SBarry Smith     levels++;
4883aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
4893aeaf226SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
4903aeaf226SBarry Smith       n    = levels;
4913aeaf226SBarry Smith       ierr = PCSetFromOptions(pc);CHKERRQ(ierr);  /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
4923aeaf226SBarry Smith       mglevels =  mg->levels;
4933aeaf226SBarry Smith     }
4943aeaf226SBarry Smith   }
4953aeaf226SBarry Smith 
496f3fbd535SBarry Smith 
497f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
498f3fbd535SBarry Smith   /* so use those from global PC */
499f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
500f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
501f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
502f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr);
5031338a6b9SJed Brown   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) {
504f3fbd535SBarry 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);
505f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
506f3fbd535SBarry Smith   }
507f3fbd535SBarry Smith 
5082d2b81a6SBarry Smith   if (pc->dm && !pc->setupcalled) {
5092d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
5102e499ae9SBarry Smith     Mat p;
51173250ac0SBarry Smith     Vec rscale;
512218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
513218a07d4SBarry Smith     dms[n-1] = pc->dm;
514218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
515218a07d4SBarry Smith       ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr);
5163c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
5173c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
51811629dbeSBarry Smith       ierr = DMSetFunction(dms[i],0);
51911629dbeSBarry Smith       ierr = DMSetInitialGuess(dms[i],0);
52024c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
521e727c939SJed Brown 	ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
5222d2b81a6SBarry Smith 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
52300ac8be1SBarry Smith 	if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
52473250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
5256bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
526218a07d4SBarry Smith       }
52724c3aa18SBarry Smith     }
5282d2b81a6SBarry Smith 
529218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
5306bf464f9SBarry Smith       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
531218a07d4SBarry Smith     }
5322d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
533cab2e9ccSBarry Smith 
534cab2e9ccSBarry Smith     /* finest smoother also gets DM but it is not active */
535cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
536cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
537218a07d4SBarry Smith   }
538218a07d4SBarry Smith 
539c91913d3SJed Brown   if (mg->galerkin == 1) {
540f3fbd535SBarry Smith     Mat B;
541f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
542f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
543f3fbd535SBarry Smith     if (!pc->setupcalled) {
544f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
545f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
546f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
547f3fbd535SBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
548f3fbd535SBarry Smith         dB   = B;
549f3fbd535SBarry Smith       }
550cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
551f3fbd535SBarry Smith     } else {
552f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
553f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
554f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
555f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
556f3fbd535SBarry Smith         dB   = B;
557f3fbd535SBarry Smith       }
558f3fbd535SBarry Smith     }
559cab2e9ccSBarry Smith   } else if (pc->dm && pc->dm->x) {
560cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
561cab2e9ccSBarry Smith     for (i=n-2;i>-1; i--) {
562*c88c5224SJed Brown       Mat R;
563*c88c5224SJed Brown       Vec rscale;
564cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
565cab2e9ccSBarry Smith         Vec *vecs;
566cab2e9ccSBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr);
567cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
568cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
569cab2e9ccSBarry Smith       }
570*c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
571*c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
572*c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
573*c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
574cab2e9ccSBarry Smith     }
575f3fbd535SBarry Smith   }
576f3fbd535SBarry Smith 
577f3fbd535SBarry Smith   if (!pc->setupcalled) {
578f3fbd535SBarry Smith     for (i=0; i<n; i++) {
579f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
580f3fbd535SBarry Smith     }
581f3fbd535SBarry Smith     for (i=1; i<n; i++) {
582f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
583f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
584f3fbd535SBarry Smith       }
585f3fbd535SBarry Smith     }
586f3fbd535SBarry Smith     for (i=1; i<n; i++) {
587*c88c5224SJed Brown       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
588*c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
589f3fbd535SBarry Smith     }
590f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
591f3fbd535SBarry Smith       if (!mglevels[i]->b) {
592f3fbd535SBarry Smith         Vec *vec;
593f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
594f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
5956bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
596f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
597f3fbd535SBarry Smith       }
598f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
599f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
600f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
6016bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
602f3fbd535SBarry Smith       }
603f3fbd535SBarry Smith       if (!mglevels[i]->x) {
604f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
605f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
6066bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
607f3fbd535SBarry Smith       }
608f3fbd535SBarry Smith     }
609f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
610f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
611f3fbd535SBarry Smith       Vec *vec;
612f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
613f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
6146bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
615f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
616f3fbd535SBarry Smith     }
617f3fbd535SBarry Smith   }
618f3fbd535SBarry Smith 
619f3fbd535SBarry Smith 
620f3fbd535SBarry Smith   for (i=1; i<n; i++) {
621f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
622f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
623f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
624f3fbd535SBarry Smith     }
62563e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
626f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
62763e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
628d42688cbSBarry Smith     if (!mglevels[i]->residual) {
629d42688cbSBarry Smith       Mat mat;
630d42688cbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
631d42688cbSBarry Smith       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
632d42688cbSBarry Smith     }
633f3fbd535SBarry Smith   }
634f3fbd535SBarry Smith   for (i=1; i<n; i++) {
635f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
636f3fbd535SBarry Smith       Mat          downmat,downpmat;
637f3fbd535SBarry Smith       MatStructure matflag;
638ace3abfcSBarry Smith       PetscBool    opsset;
639f3fbd535SBarry Smith 
640f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
641f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
642f3fbd535SBarry Smith       if (!opsset) {
643f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
644f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
645f3fbd535SBarry Smith       }
646f3fbd535SBarry Smith 
647f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
64863e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
649f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
65063e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
651f3fbd535SBarry Smith     }
652f3fbd535SBarry Smith   }
653f3fbd535SBarry Smith 
654f3fbd535SBarry Smith   /*
655f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
656f3fbd535SBarry Smith   */
657f3fbd535SBarry Smith   ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
658f3fbd535SBarry Smith   if (preonly) {
659f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
660f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
661f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
662fd1303e9SJungho Lee     ierr = PetscTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
663fd1303e9SJungho Lee     if (!lu && !redundant && !cholesky && !svd) {
664f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
665f3fbd535SBarry Smith     }
666f3fbd535SBarry Smith   }
667f3fbd535SBarry Smith 
668f3fbd535SBarry Smith   if (!pc->setupcalled) {
669f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
670f3fbd535SBarry Smith   }
671f3fbd535SBarry Smith 
67263e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
673f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
67463e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
675f3fbd535SBarry Smith 
676f3fbd535SBarry Smith   /*
677f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
678e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
679f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
680f3fbd535SBarry Smith 
681f3fbd535SBarry Smith    Only support one or the other at the same time.
682f3fbd535SBarry Smith   */
683f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
684acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
685f3fbd535SBarry Smith   if (dump) {
686f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
687f3fbd535SBarry Smith   }
688f3fbd535SBarry Smith   dump = PETSC_FALSE;
689f3fbd535SBarry Smith #endif
690acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
691f3fbd535SBarry Smith   if (dump) {
692f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
693f3fbd535SBarry Smith   }
694f3fbd535SBarry Smith 
695f3fbd535SBarry Smith   if (viewer) {
696f3fbd535SBarry Smith     for (i=1; i<n; i++) {
697f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
698f3fbd535SBarry Smith     }
699f3fbd535SBarry Smith     for (i=0; i<n; i++) {
700f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
701f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
702f3fbd535SBarry Smith     }
703f3fbd535SBarry Smith   }
704f3fbd535SBarry Smith   PetscFunctionReturn(0);
705f3fbd535SBarry Smith }
706f3fbd535SBarry Smith 
707f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
708f3fbd535SBarry Smith 
709f3fbd535SBarry Smith #undef __FUNCT__
7109dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
7114b9ad928SBarry Smith /*@
71297177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
7134b9ad928SBarry Smith 
7144b9ad928SBarry Smith    Not Collective
7154b9ad928SBarry Smith 
7164b9ad928SBarry Smith    Input Parameter:
7174b9ad928SBarry Smith .  pc - the preconditioner context
7184b9ad928SBarry Smith 
7194b9ad928SBarry Smith    Output parameter:
7204b9ad928SBarry Smith .  levels - the number of levels
7214b9ad928SBarry Smith 
7224b9ad928SBarry Smith    Level: advanced
7234b9ad928SBarry Smith 
7244b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
7254b9ad928SBarry Smith 
72697177400SBarry Smith .seealso: PCMGSetLevels()
7274b9ad928SBarry Smith @*/
7287087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
7294b9ad928SBarry Smith {
730f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
7314b9ad928SBarry Smith 
7324b9ad928SBarry Smith   PetscFunctionBegin;
7330700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
7344482741eSBarry Smith   PetscValidIntPointer(levels,2);
735f3fbd535SBarry Smith   *levels = mg->nlevels;
7364b9ad928SBarry Smith   PetscFunctionReturn(0);
7374b9ad928SBarry Smith }
7384b9ad928SBarry Smith 
7394b9ad928SBarry Smith #undef __FUNCT__
7409dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
7414b9ad928SBarry Smith /*@
74297177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
7434b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
7444b9ad928SBarry Smith 
745ad4df100SBarry Smith    Logically Collective on PC
7464b9ad928SBarry Smith 
7474b9ad928SBarry Smith    Input Parameters:
7484b9ad928SBarry Smith +  pc - the preconditioner context
7499dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
7509dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
7514b9ad928SBarry Smith 
7524b9ad928SBarry Smith    Options Database Key:
7534b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
7544b9ad928SBarry Smith    additive, full, kaskade
7554b9ad928SBarry Smith 
7564b9ad928SBarry Smith    Level: advanced
7574b9ad928SBarry Smith 
7584b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
7594b9ad928SBarry Smith 
76097177400SBarry Smith .seealso: PCMGSetLevels()
7614b9ad928SBarry Smith @*/
7627087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
7634b9ad928SBarry Smith {
764f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
7654b9ad928SBarry Smith 
7664b9ad928SBarry Smith   PetscFunctionBegin;
7670700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
768c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
76931567311SBarry Smith   mg->am = form;
7709dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
7714b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
7724b9ad928SBarry Smith   PetscFunctionReturn(0);
7734b9ad928SBarry Smith }
7744b9ad928SBarry Smith 
7754b9ad928SBarry Smith #undef __FUNCT__
7760d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
7774b9ad928SBarry Smith /*@
7780d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
7794b9ad928SBarry Smith    complicated cycling.
7804b9ad928SBarry Smith 
781ad4df100SBarry Smith    Logically Collective on PC
7824b9ad928SBarry Smith 
7834b9ad928SBarry Smith    Input Parameters:
784c2be2410SBarry Smith +  pc - the multigrid context
7850d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
7864b9ad928SBarry Smith 
7874b9ad928SBarry Smith    Options Database Key:
7880d353602SBarry Smith $  -pc_mg_cycle_type v or w
7894b9ad928SBarry Smith 
7904b9ad928SBarry Smith    Level: advanced
7914b9ad928SBarry Smith 
7924b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7934b9ad928SBarry Smith 
7940d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
7954b9ad928SBarry Smith @*/
7967087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
7974b9ad928SBarry Smith {
798f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
799f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
80079416396SBarry Smith   PetscInt     i,levels;
8014b9ad928SBarry Smith 
8024b9ad928SBarry Smith   PetscFunctionBegin;
8030700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
804e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
805c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
806f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8074b9ad928SBarry Smith 
8084b9ad928SBarry Smith   for (i=0; i<levels; i++) {
809f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
8104b9ad928SBarry Smith   }
8114b9ad928SBarry Smith   PetscFunctionReturn(0);
8124b9ad928SBarry Smith }
8134b9ad928SBarry Smith 
8144b9ad928SBarry Smith #undef __FUNCT__
8158cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
8168cc2d5dfSBarry Smith /*@
8178cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
8188cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
8198cc2d5dfSBarry Smith 
820ad4df100SBarry Smith    Logically Collective on PC
8218cc2d5dfSBarry Smith 
8228cc2d5dfSBarry Smith    Input Parameters:
8238cc2d5dfSBarry Smith +  pc - the multigrid context
8248cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
8258cc2d5dfSBarry Smith 
8268cc2d5dfSBarry Smith    Options Database Key:
8278cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
8288cc2d5dfSBarry Smith 
8298cc2d5dfSBarry Smith    Level: advanced
8308cc2d5dfSBarry Smith 
8318cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
8328cc2d5dfSBarry Smith 
8338cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8348cc2d5dfSBarry Smith 
8358cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
8368cc2d5dfSBarry Smith @*/
8377087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
8388cc2d5dfSBarry Smith {
839f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
840f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
8418cc2d5dfSBarry Smith   PetscInt     i,levels;
8428cc2d5dfSBarry Smith 
8438cc2d5dfSBarry Smith   PetscFunctionBegin;
8440700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
845e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
846c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
847f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8488cc2d5dfSBarry Smith 
8498cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
85031567311SBarry Smith     mg->cyclesperpcapply  = n;
8518cc2d5dfSBarry Smith   }
8528cc2d5dfSBarry Smith   PetscFunctionReturn(0);
8538cc2d5dfSBarry Smith }
8548cc2d5dfSBarry Smith 
8558cc2d5dfSBarry Smith #undef __FUNCT__
8569dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
857c2be2410SBarry Smith /*@
85897177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
859c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
860c2be2410SBarry Smith 
861ad4df100SBarry Smith    Logically Collective on PC
862c2be2410SBarry Smith 
863c2be2410SBarry Smith    Input Parameters:
864c91913d3SJed Brown +  pc - the multigrid context
865c91913d3SJed Brown -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
866c2be2410SBarry Smith 
867c2be2410SBarry Smith    Options Database Key:
868c2be2410SBarry Smith $  -pc_mg_galerkin
869c2be2410SBarry Smith 
870c2be2410SBarry Smith    Level: intermediate
871c2be2410SBarry Smith 
872c2be2410SBarry Smith .keywords: MG, set, Galerkin
873c2be2410SBarry Smith 
8743fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
8753fc8bf9cSBarry Smith 
876c2be2410SBarry Smith @*/
877c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
878c2be2410SBarry Smith {
879f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
880c2be2410SBarry Smith 
881c2be2410SBarry Smith   PetscFunctionBegin;
8820700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
883c91913d3SJed Brown   mg->galerkin = (PetscInt)use;
884c2be2410SBarry Smith   PetscFunctionReturn(0);
885c2be2410SBarry Smith }
886c2be2410SBarry Smith 
887c2be2410SBarry Smith #undef __FUNCT__
8883fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
8893fc8bf9cSBarry Smith /*@
8903fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
8913fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
8923fc8bf9cSBarry Smith 
8933fc8bf9cSBarry Smith    Not Collective
8943fc8bf9cSBarry Smith 
8953fc8bf9cSBarry Smith    Input Parameter:
8963fc8bf9cSBarry Smith .  pc - the multigrid context
8973fc8bf9cSBarry Smith 
8983fc8bf9cSBarry Smith    Output Parameter:
8993fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
9003fc8bf9cSBarry Smith 
9013fc8bf9cSBarry Smith    Options Database Key:
9023fc8bf9cSBarry Smith $  -pc_mg_galerkin
9033fc8bf9cSBarry Smith 
9043fc8bf9cSBarry Smith    Level: intermediate
9053fc8bf9cSBarry Smith 
9063fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
9073fc8bf9cSBarry Smith 
9083fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
9093fc8bf9cSBarry Smith 
9103fc8bf9cSBarry Smith @*/
9117087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
9123fc8bf9cSBarry Smith {
913f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
9143fc8bf9cSBarry Smith 
9153fc8bf9cSBarry Smith   PetscFunctionBegin;
9160700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
91713fdb3acSJose Roman   *galerkin = (PetscBool)mg->galerkin;
9183fc8bf9cSBarry Smith   PetscFunctionReturn(0);
9193fc8bf9cSBarry Smith }
9203fc8bf9cSBarry Smith 
9213fc8bf9cSBarry Smith #undef __FUNCT__
9229dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
9234b9ad928SBarry Smith /*@
92497177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
92597177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
9264b9ad928SBarry Smith    pre-smoothing steps on different levels.
9274b9ad928SBarry Smith 
928ad4df100SBarry Smith    Logically Collective on PC
9294b9ad928SBarry Smith 
9304b9ad928SBarry Smith    Input Parameters:
9314b9ad928SBarry Smith +  mg - the multigrid context
9324b9ad928SBarry Smith -  n - the number of smoothing steps
9334b9ad928SBarry Smith 
9344b9ad928SBarry Smith    Options Database Key:
9354b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
9364b9ad928SBarry Smith 
9374b9ad928SBarry Smith    Level: advanced
9384b9ad928SBarry Smith 
9394b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
9404b9ad928SBarry Smith 
94197177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
9424b9ad928SBarry Smith @*/
9437087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
9444b9ad928SBarry Smith {
945f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
946f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9476849ba73SBarry Smith   PetscErrorCode ierr;
94879416396SBarry Smith   PetscInt       i,levels;
9494b9ad928SBarry Smith 
9504b9ad928SBarry Smith   PetscFunctionBegin;
9510700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
952e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
953c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
954f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9554b9ad928SBarry Smith 
956b05257ddSBarry Smith   for (i=1; i<levels; i++) {
9574b9ad928SBarry Smith     /* make sure smoother up and down are different */
95897177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
959f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
96031567311SBarry Smith     mg->default_smoothd = n;
9614b9ad928SBarry Smith   }
9624b9ad928SBarry Smith   PetscFunctionReturn(0);
9634b9ad928SBarry Smith }
9644b9ad928SBarry Smith 
9654b9ad928SBarry Smith #undef __FUNCT__
9669dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
9674b9ad928SBarry Smith /*@
96897177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
96997177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
9704b9ad928SBarry Smith    post-smoothing steps on different levels.
9714b9ad928SBarry Smith 
972ad4df100SBarry Smith    Logically Collective on PC
9734b9ad928SBarry Smith 
9744b9ad928SBarry Smith    Input Parameters:
9754b9ad928SBarry Smith +  mg - the multigrid context
9764b9ad928SBarry Smith -  n - the number of smoothing steps
9774b9ad928SBarry Smith 
9784b9ad928SBarry Smith    Options Database Key:
9794b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
9804b9ad928SBarry Smith 
9814b9ad928SBarry Smith    Level: advanced
9824b9ad928SBarry Smith 
9834b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
984a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
9854b9ad928SBarry Smith 
9864b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
9874b9ad928SBarry Smith 
98897177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
9894b9ad928SBarry Smith @*/
9907087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
9914b9ad928SBarry Smith {
992f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
993f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9946849ba73SBarry Smith   PetscErrorCode ierr;
99579416396SBarry Smith   PetscInt       i,levels;
9964b9ad928SBarry Smith 
9974b9ad928SBarry Smith   PetscFunctionBegin;
9980700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
999e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1000c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1001f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10024b9ad928SBarry Smith 
10034b9ad928SBarry Smith   for (i=1; i<levels; i++) {
10044b9ad928SBarry Smith     /* make sure smoother up and down are different */
100597177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1006f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
100731567311SBarry Smith     mg->default_smoothu = n;
10084b9ad928SBarry Smith   }
10094b9ad928SBarry Smith   PetscFunctionReturn(0);
10104b9ad928SBarry Smith }
10114b9ad928SBarry Smith 
10124b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
10134b9ad928SBarry Smith 
10143b09bd56SBarry Smith /*MC
1015ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
10163b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
10173b09bd56SBarry Smith 
10183b09bd56SBarry Smith    Options Database Keys:
10193b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
10200d353602SBarry Smith .  -pc_mg_cycles v or w
102179416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
10223b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
10233b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
10243b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
10253b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
102668eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
10273b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1028e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
10293b09bd56SBarry Smith 
103024c3aa18SBarry 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
10313b09bd56SBarry Smith 
10323b09bd56SBarry Smith    Level: intermediate
10333b09bd56SBarry Smith 
10348f87f92bSBarry Smith    Concepts: multigrid/multilevel
10353b09bd56SBarry Smith 
103624c3aa18SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC
10370d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
103897177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
103997177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
10400d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
10413b09bd56SBarry Smith M*/
10423b09bd56SBarry Smith 
10434b9ad928SBarry Smith EXTERN_C_BEGIN
10444b9ad928SBarry Smith #undef __FUNCT__
10454b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
10467087cfbeSBarry Smith PetscErrorCode  PCCreate_MG(PC pc)
10474b9ad928SBarry Smith {
1048f3fbd535SBarry Smith   PC_MG          *mg;
1049f3fbd535SBarry Smith   PetscErrorCode ierr;
1050f3fbd535SBarry Smith 
10514b9ad928SBarry Smith   PetscFunctionBegin;
1052f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1053f3fbd535SBarry Smith   pc->data    = (void*)mg;
1054f3fbd535SBarry Smith   mg->nlevels = -1;
1055f3fbd535SBarry Smith 
10564b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
10574b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1058a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
10594b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
10604b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
10614b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
10624b9ad928SBarry Smith   PetscFunctionReturn(0);
10634b9ad928SBarry Smith }
10644b9ad928SBarry Smith EXTERN_C_END
1065