xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 1cd381d66697d5ca2812fcd8e67c7fc508f6cf09)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscpcmg.h" I*/
64b9ad928SBarry Smith 
74b9ad928SBarry Smith 
84b9ad928SBarry Smith #undef __FUNCT__
99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private"
1031567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
114b9ad928SBarry Smith {
1231567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1331567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
146849ba73SBarry Smith   PetscErrorCode ierr;
1531567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
184b9ad928SBarry Smith 
1963e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2031567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
2163e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2231567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2363e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2431567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2563e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
264b9ad928SBarry Smith 
274b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
2831567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
294b9ad928SBarry Smith       PetscReal rnorm;
3031567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
314b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3270441072SBarry Smith         if (rnorm < mg->abstol) {
334d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
341e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
354b9ad928SBarry Smith         } else {
364d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
371e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr);
384b9ad928SBarry Smith         }
394b9ad928SBarry Smith         PetscFunctionReturn(0);
404b9ad928SBarry Smith       }
414b9ad928SBarry Smith     }
424b9ad928SBarry Smith 
4331567311SBarry Smith     mgc = *(mglevelsin - 1);
4463e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4531567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
4663e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
47efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
484b9ad928SBarry Smith     while (cycles--) {
4931567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
504b9ad928SBarry Smith     }
5163e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5231567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5363e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5463e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5531567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
5663e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith   }
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
614b9ad928SBarry Smith #undef __FUNCT__
624b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
63ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool  zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
644b9ad928SBarry Smith {
65f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
66f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
67dfbe8321SBarry Smith   PetscErrorCode ierr;
68f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
694b9ad928SBarry Smith 
704b9ad928SBarry Smith   PetscFunctionBegin;
71f3fbd535SBarry Smith   mglevels[levels-1]->b    = b;
72f3fbd535SBarry Smith   mglevels[levels-1]->x    = x;
734b9ad928SBarry Smith 
7431567311SBarry Smith   mg->rtol = rtol;
7531567311SBarry Smith   mg->abstol = abstol;
7631567311SBarry Smith   mg->dtol = dtol;
774b9ad928SBarry Smith   if (rtol) {
784b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
794b9ad928SBarry Smith     PetscReal rnorm;
807319c654SBarry Smith     if (zeroguess) {
817319c654SBarry Smith       ierr               = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
827319c654SBarry Smith     } else {
83f3fbd535SBarry Smith       ierr               = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
844b9ad928SBarry Smith       ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
857319c654SBarry Smith     }
8631567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
8770441072SBarry Smith   } else if (abstol) {
8831567311SBarry Smith     mg->ttol = abstol;
894b9ad928SBarry Smith   } else {
9031567311SBarry Smith     mg->ttol = 0.0;
914b9ad928SBarry Smith   }
924b9ad928SBarry Smith 
934d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
94336babb1SJed Brown      stop prematurely due to small residual */
954d0a8057SBarry Smith   for (i=1; i<levels; i++) {
96f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
97f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
98f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
994b9ad928SBarry Smith     }
1004d0a8057SBarry Smith   }
1014d0a8057SBarry Smith 
1024d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1034d0a8057SBarry Smith   for (i=0; i<its; i++) {
104f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1054d0a8057SBarry Smith     if (*reason) break;
1064d0a8057SBarry Smith   }
1074d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1084d0a8057SBarry Smith   *outits = i;
1094b9ad928SBarry Smith   PetscFunctionReturn(0);
1104b9ad928SBarry Smith }
1114b9ad928SBarry Smith 
1124b9ad928SBarry Smith #undef __FUNCT__
1133aeaf226SBarry Smith #define __FUNCT__ "PCReset_MG"
1143aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1153aeaf226SBarry Smith {
1163aeaf226SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1173aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1183aeaf226SBarry Smith   PetscErrorCode ierr;
1193aeaf226SBarry Smith   PetscInt       i,n;
1203aeaf226SBarry Smith 
1213aeaf226SBarry Smith   PetscFunctionBegin;
1223aeaf226SBarry Smith   if (mglevels) {
1233aeaf226SBarry Smith     n = mglevels[0]->levels;
1243aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1253aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1263aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1273aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
1283aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1293aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
13073250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1313aeaf226SBarry Smith     }
1323aeaf226SBarry Smith 
1333aeaf226SBarry Smith     for (i=0; i<n; i++) {
1343aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1353aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1363aeaf226SBarry Smith 	ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1373aeaf226SBarry Smith       }
1383aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1393aeaf226SBarry Smith     }
1403aeaf226SBarry Smith   }
1413aeaf226SBarry Smith   PetscFunctionReturn(0);
1423aeaf226SBarry Smith }
1433aeaf226SBarry Smith 
1443aeaf226SBarry Smith #undef __FUNCT__
1459dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
1464b9ad928SBarry Smith /*@C
14797177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1484b9ad928SBarry Smith    Must be called before any other MG routine.
1494b9ad928SBarry Smith 
150ad4df100SBarry Smith    Logically Collective on PC
1514b9ad928SBarry Smith 
1524b9ad928SBarry Smith    Input Parameters:
1534b9ad928SBarry Smith +  pc - the preconditioner context
1544b9ad928SBarry Smith .  levels - the number of levels
1554b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1564b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
1574b9ad928SBarry Smith 
1584b9ad928SBarry Smith    Level: intermediate
1594b9ad928SBarry Smith 
1604b9ad928SBarry Smith    Notes:
1614b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1624b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1634b9ad928SBarry Smith 
1644b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1654b9ad928SBarry Smith 
16697177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1674b9ad928SBarry Smith @*/
1687087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1694b9ad928SBarry Smith {
170dfbe8321SBarry Smith   PetscErrorCode ierr;
171f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
172f3fbd535SBarry Smith   MPI_Comm       comm = ((PetscObject)pc)->comm;
1733aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
174f3fbd535SBarry Smith   PetscInt       i;
175f3fbd535SBarry Smith   PetscMPIInt    size;
176f3fbd535SBarry Smith   const char     *prefix;
177f3fbd535SBarry Smith   PC             ipc;
1783aeaf226SBarry Smith   PetscInt       n;
1794b9ad928SBarry Smith 
1804b9ad928SBarry Smith   PetscFunctionBegin;
1810700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
182548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
183548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1843aeaf226SBarry Smith   if (mglevels) {
1853aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
1863aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
1873aeaf226SBarry Smith     n = mglevels[0]->levels;
1883aeaf226SBarry Smith     for (i=0; i<n; i++) {
1893aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1903aeaf226SBarry Smith 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
1913aeaf226SBarry Smith       }
1923aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
1933aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
1943aeaf226SBarry Smith     }
1953aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
1963aeaf226SBarry Smith   }
197f3fbd535SBarry Smith 
198f3fbd535SBarry Smith   mg->nlevels      = levels;
199f3fbd535SBarry Smith 
200f3fbd535SBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr);
201f3fbd535SBarry Smith   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
202f3fbd535SBarry Smith 
203f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
204f3fbd535SBarry Smith 
205a9db87a2SMatthew G Knepley   mg->stageApply = 0;
206f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
207f3fbd535SBarry Smith     ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr);
208f3fbd535SBarry Smith     mglevels[i]->level           = i;
209f3fbd535SBarry Smith     mglevels[i]->levels          = levels;
210f3fbd535SBarry Smith     mglevels[i]->cycles          = PC_MG_CYCLE_V;
211336babb1SJed Brown     mg->default_smoothu = 2;
212336babb1SJed Brown     mg->default_smoothd = 2;
21363e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
21463e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
21563e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
21663e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
217f3fbd535SBarry Smith 
218f3fbd535SBarry Smith     if (comms) comm = comms[i];
219f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
220336babb1SJed Brown     ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
221336babb1SJed Brown     ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
222336babb1SJed Brown     ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
223336babb1SJed Brown     ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
224336babb1SJed Brown     ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
225f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
226336babb1SJed Brown     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i?mg->default_smoothd:1);CHKERRQ(ierr);
227f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
228f3fbd535SBarry Smith 
229f3fbd535SBarry Smith     /* do special stuff for coarse grid */
230f3fbd535SBarry Smith     if (!i && levels > 1) {
231f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
232f3fbd535SBarry Smith 
233f3fbd535SBarry Smith       /* coarse solve is (redundant) LU by default */
234f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
235f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
236f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
237f3fbd535SBarry Smith       if (size > 1) {
238f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
239f3fbd535SBarry Smith       } else {
240f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
241f3fbd535SBarry Smith       }
242f3fbd535SBarry Smith 
243f3fbd535SBarry Smith     } else {
244f3fbd535SBarry Smith       char tprefix[128];
245f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
246f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
247f3fbd535SBarry Smith     }
248f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr);
249f3fbd535SBarry Smith     mglevels[i]->smoothu    = mglevels[i]->smoothd;
25031567311SBarry Smith     mg->rtol                = 0.0;
25131567311SBarry Smith     mg->abstol              = 0.0;
25231567311SBarry Smith     mg->dtol                = 0.0;
25331567311SBarry Smith     mg->ttol                = 0.0;
25431567311SBarry Smith     mg->cyclesperpcapply    = 1;
255f3fbd535SBarry Smith   }
25631567311SBarry Smith   mg->am          = PC_MG_MULTIPLICATIVE;
257f3fbd535SBarry Smith   mg->levels      = mglevels;
2584b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
2594b9ad928SBarry Smith   PetscFunctionReturn(0);
2604b9ad928SBarry Smith }
2614b9ad928SBarry Smith 
262c07bf074SBarry Smith 
263c07bf074SBarry Smith #undef __FUNCT__
264c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
265c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
266c07bf074SBarry Smith {
267c07bf074SBarry Smith   PetscErrorCode ierr;
268a06653b4SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
269a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
270a06653b4SBarry Smith   PetscInt       i,n;
271c07bf074SBarry Smith 
272c07bf074SBarry Smith   PetscFunctionBegin;
273a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
274a06653b4SBarry Smith   if (mglevels) {
275a06653b4SBarry Smith     n = mglevels[0]->levels;
276a06653b4SBarry Smith     for (i=0; i<n; i++) {
277a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2786bf464f9SBarry Smith 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
279a06653b4SBarry Smith       }
2806bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
281a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
282a06653b4SBarry Smith     }
283a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
284a06653b4SBarry Smith   }
285c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
286f3fbd535SBarry Smith   PetscFunctionReturn(0);
287f3fbd535SBarry Smith }
288f3fbd535SBarry Smith 
289f3fbd535SBarry Smith 
290f3fbd535SBarry Smith 
29109573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
29209573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
29309573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
294f3fbd535SBarry Smith 
295f3fbd535SBarry Smith /*
296f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
297f3fbd535SBarry Smith              or full cycle of multigrid.
298f3fbd535SBarry Smith 
299f3fbd535SBarry Smith   Note:
300f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
301f3fbd535SBarry Smith */
302f3fbd535SBarry Smith #undef __FUNCT__
303f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
304f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
305f3fbd535SBarry Smith {
306f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
307f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
308f3fbd535SBarry Smith   PetscErrorCode ierr;
309f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
310f3fbd535SBarry Smith 
311f3fbd535SBarry Smith   PetscFunctionBegin;
312a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
313e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
314298cc208SBarry Smith   for (i=0; i<levels; i++) {
315e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
316e1d8e5deSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
317298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
318e1d8e5deSBarry Smith     }
319e1d8e5deSBarry Smith   }
320e1d8e5deSBarry Smith 
321f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
322f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
32331567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
324f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
32531567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
326f3fbd535SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
327f3fbd535SBarry Smith     }
328f3fbd535SBarry Smith   }
32931567311SBarry Smith   else if (mg->am == PC_MG_ADDITIVE) {
33031567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
331f3fbd535SBarry Smith   }
33231567311SBarry Smith   else if (mg->am == PC_MG_KASKADE) {
33331567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
334f3fbd535SBarry Smith   }
335f3fbd535SBarry Smith   else {
336f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
337f3fbd535SBarry Smith   }
338a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
339f3fbd535SBarry Smith   PetscFunctionReturn(0);
340f3fbd535SBarry Smith }
341f3fbd535SBarry Smith 
342f3fbd535SBarry Smith 
343f3fbd535SBarry Smith #undef __FUNCT__
344f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
345f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc)
346f3fbd535SBarry Smith {
347f3fbd535SBarry Smith   PetscErrorCode ierr;
348f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
349c91913d3SJed Brown   PetscBool      flg,set;
350f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
351f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
352f3fbd535SBarry Smith   PCMGType       mgtype;
353f3fbd535SBarry Smith   PCMGCycleType  mgctype;
354f3fbd535SBarry Smith 
355f3fbd535SBarry Smith   PetscFunctionBegin;
356f3fbd535SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
3573aeaf226SBarry Smith     if (!mg->levels) {
358f3fbd535SBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
359eb3f98d2SBarry Smith       if (!flg && pc->dm) {
360eb3f98d2SBarry Smith         ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
361eb3f98d2SBarry Smith         levels++;
3623aeaf226SBarry Smith         mg->usedmfornumberoflevels = PETSC_TRUE;
363eb3f98d2SBarry Smith       }
364f3fbd535SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
365f3fbd535SBarry Smith     }
3663aeaf226SBarry Smith     mglevels = mg->levels;
3673aeaf226SBarry Smith 
368f3fbd535SBarry Smith     mgctype = (PCMGCycleType) mglevels[0]->cycles;
369f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
370f3fbd535SBarry Smith     if (flg) {
371f3fbd535SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
372f3fbd535SBarry Smith     };
373f3fbd535SBarry Smith     flg  = PETSC_FALSE;
374c91913d3SJed Brown     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr);
375c91913d3SJed Brown     if (set) {
376c91913d3SJed Brown       ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr);
377f3fbd535SBarry Smith     }
378f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
379f3fbd535SBarry Smith     if (flg) {
380f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
381f3fbd535SBarry Smith     }
382f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
383f3fbd535SBarry Smith     if (flg) {
384f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
385f3fbd535SBarry Smith     }
38631567311SBarry Smith     mgtype = mg->am;
387f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
388f3fbd535SBarry Smith     if (flg) {
389f3fbd535SBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
390f3fbd535SBarry Smith     }
39131567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
39231567311SBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
393f3fbd535SBarry Smith       if (flg) {
394f3fbd535SBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
395f3fbd535SBarry Smith       }
396f3fbd535SBarry Smith     }
397f3fbd535SBarry Smith     flg  = PETSC_FALSE;
398acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
399f3fbd535SBarry Smith     if (flg) {
400f3fbd535SBarry Smith       PetscInt i;
401f3fbd535SBarry Smith       char     eventname[128];
40263e6d426SJed Brown       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
403f3fbd535SBarry Smith       levels = mglevels[0]->levels;
404f3fbd535SBarry Smith       for (i=0; i<levels; i++) {
405f3fbd535SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
40663e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
407f3fbd535SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
40863e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
409f3fbd535SBarry Smith         if (i) {
410f3fbd535SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
41163e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
412f3fbd535SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
41363e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
414f3fbd535SBarry Smith         }
415f3fbd535SBarry Smith       }
416eec35531SMatthew G Knepley 
417ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
418eec35531SMatthew G Knepley       {
419eec35531SMatthew G Knepley         const char   *sname = "MG Apply";
420eec35531SMatthew G Knepley         PetscStageLog stageLog;
421eec35531SMatthew G Knepley         PetscInt      st;
422eec35531SMatthew G Knepley 
423eec35531SMatthew G Knepley         PetscFunctionBegin;
424eec35531SMatthew G Knepley         ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
425eec35531SMatthew G Knepley         for (st = 0; st < stageLog->numStages; ++st) {
426eec35531SMatthew G Knepley           PetscBool same;
427eec35531SMatthew G Knepley 
428eec35531SMatthew G Knepley           ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
429eec35531SMatthew G Knepley           if (same) {mg->stageApply = st;}
430eec35531SMatthew G Knepley         }
431eec35531SMatthew G Knepley         if (!mg->stageApply) {
432eec35531SMatthew G Knepley           ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
433eec35531SMatthew G Knepley         }
434eec35531SMatthew G Knepley       }
435ec60ca73SSatish Balay #endif
436f3fbd535SBarry Smith     }
437f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
438f3fbd535SBarry Smith   PetscFunctionReturn(0);
439f3fbd535SBarry Smith }
440f3fbd535SBarry Smith 
4416a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4426a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
443f3fbd535SBarry Smith 
444f3fbd535SBarry Smith #undef __FUNCT__
445f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
446f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
447f3fbd535SBarry Smith {
448f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
449f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
450f3fbd535SBarry Smith   PetscErrorCode ierr;
451f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
452219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
453f3fbd535SBarry Smith 
454f3fbd535SBarry Smith   PetscFunctionBegin;
455251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4565b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
457219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
458f3fbd535SBarry Smith   if (iascii) {
45931567311SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,(mglevels[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr);
46031567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
46131567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
462f3fbd535SBarry Smith     }
463218a07d4SBarry Smith     if (mg->galerkin) {
464f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4654f66f45eSBarry Smith     } else {
4664f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
467f3fbd535SBarry Smith     }
468f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
469f3fbd535SBarry Smith       if (!i) {
47089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
471f3fbd535SBarry Smith       } else {
47289cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
473f3fbd535SBarry Smith       }
474f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
475f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
476f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
477f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
478f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
479f3fbd535SBarry Smith       } else if (i){
48089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
481f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
482f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
483f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
484f3fbd535SBarry Smith       }
485f3fbd535SBarry Smith     }
4865b0b0462SBarry Smith   } else if (isbinary) {
4875b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
4885b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
4895b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
4905b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
4915b0b0462SBarry Smith       }
4925b0b0462SBarry Smith     }
493219b1045SBarry Smith   } else if (isdraw) {
494219b1045SBarry Smith     PetscDraw draw;
495b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
496219b1045SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
497219b1045SBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
498219b1045SBarry Smith     ierr = PetscDrawStringGetSize(draw,PETSC_NULL,&th);CHKERRQ(ierr);
499219b1045SBarry Smith     bottom = y - th;
500219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
501b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
502219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
503219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
504219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
505b4375e8dSPeter Brune       } else {
506b4375e8dSPeter Brune         w = 0.5*PetscMin(1.0-x,x);
507b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
508b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
509b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
510b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
511b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
512b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
513b4375e8dSPeter Brune       }
514*1cd381d6SBarry Smith       ierr = PetscDrawGetBoundingBox(draw,PETSC_NULL,&bottom,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
515*1cd381d6SBarry Smith       bottom -= th;
516219b1045SBarry Smith     }
5175b0b0462SBarry Smith   }
518f3fbd535SBarry Smith   PetscFunctionReturn(0);
519f3fbd535SBarry Smith }
520f3fbd535SBarry Smith 
521b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>
522b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
523cab2e9ccSBarry Smith 
524f3fbd535SBarry Smith /*
525f3fbd535SBarry Smith     Calls setup for the KSP on each level
526f3fbd535SBarry Smith */
527f3fbd535SBarry Smith #undef __FUNCT__
528f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
529f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
530f3fbd535SBarry Smith {
531f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
532f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
533f3fbd535SBarry Smith   PetscErrorCode          ierr;
534f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
53598e04984SBarry Smith   PC                      cpc;
536649052a6SBarry Smith   PetscBool               preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset;
537f3fbd535SBarry Smith   Mat                     dA,dB;
538f3fbd535SBarry Smith   MatStructure            uflag;
539f3fbd535SBarry Smith   Vec                     tvec;
540218a07d4SBarry Smith   DM                      *dms;
541649052a6SBarry Smith   PetscViewer             viewer = 0;
542f3fbd535SBarry Smith 
543f3fbd535SBarry Smith   PetscFunctionBegin;
54401bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5453aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5463aeaf226SBarry Smith     PetscInt levels;
5473aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5483aeaf226SBarry Smith     levels++;
5493aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5503aeaf226SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
5513aeaf226SBarry Smith       n    = levels;
5523aeaf226SBarry Smith       ierr = PCSetFromOptions(pc);CHKERRQ(ierr);  /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5533aeaf226SBarry Smith       mglevels =  mg->levels;
5543aeaf226SBarry Smith     }
5553aeaf226SBarry Smith   }
55698e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5573aeaf226SBarry Smith 
558f3fbd535SBarry Smith 
559f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
560f3fbd535SBarry Smith   /* so use those from global PC */
561f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
562f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
56398e04984SBarry Smith   if (opsset) {
56498e04984SBarry Smith     Mat mmat;
56598e04984SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,PETSC_NULL,&mmat,PETSC_NULL);CHKERRQ(ierr);
56698e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
56798e04984SBarry Smith   }
56898e04984SBarry Smith   if (!opsset) {
569f3fbd535SBarry 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);
570f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
571f3fbd535SBarry Smith   }
572f3fbd535SBarry Smith 
57301bc414fSDmitry Karpeev   /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */
574ce4cda84SJed Brown   if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
5752d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
5762e499ae9SBarry Smith     Mat p;
57773250ac0SBarry Smith     Vec rscale;
578218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
579218a07d4SBarry Smith     dms[n-1] = pc->dm;
580218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
581942e3340SBarry Smith       DMKSP kdm;
5822ee06e3bSJed Brown       ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);
5833c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
5843c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
585942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
586d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
587d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
5885358d0d4SBarry Smith       kdm->ops->computerhs = PETSC_NULL;
589fe86f630SJed Brown       kdm->rhsctx          = PETSC_NULL;
59024c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
591e727c939SJed Brown 	ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
5922d2b81a6SBarry Smith 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
59300ac8be1SBarry Smith 	if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
59473250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
5956bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
596218a07d4SBarry Smith       }
59724c3aa18SBarry Smith     }
5982d2b81a6SBarry Smith 
599218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
6006bf464f9SBarry Smith       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
601218a07d4SBarry Smith     }
6022d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
603ce4cda84SJed Brown   }
604cab2e9ccSBarry Smith 
605ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
606ce4cda84SJed Brown     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
607cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
608cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
609218a07d4SBarry Smith   }
610218a07d4SBarry Smith 
611c91913d3SJed Brown   if (mg->galerkin == 1) {
612f3fbd535SBarry Smith     Mat B;
613f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
614f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
615f3fbd535SBarry Smith     if (!pc->setupcalled) {
616f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
617f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
618f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
619f3fbd535SBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
620f3fbd535SBarry Smith         dB   = B;
621f3fbd535SBarry Smith       }
622cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
623f3fbd535SBarry Smith     } else {
624f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
625f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
626f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
627f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
628f3fbd535SBarry Smith         dB   = B;
629f3fbd535SBarry Smith       }
630f3fbd535SBarry Smith     }
631ce4cda84SJed Brown   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
632cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
633cab2e9ccSBarry Smith     for (i=n-2;i>-1; i--) {
634c88c5224SJed Brown       Mat R;
635c88c5224SJed Brown       Vec rscale;
636cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
637cab2e9ccSBarry Smith         Vec *vecs;
638cab2e9ccSBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr);
639cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
640cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
641cab2e9ccSBarry Smith       }
642c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
643c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
644c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
645c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
646cab2e9ccSBarry Smith     }
647f3fbd535SBarry Smith   }
648ccceb50cSJed Brown   if (!mg->galerkin && pc->dm) {
649caa4e7f2SJed Brown     for (i=n-2;i>=0; i--) {
650ccceb50cSJed Brown       DM dmfine,dmcoarse;
651caa4e7f2SJed Brown       Mat Restrict,Inject;
652caa4e7f2SJed Brown       Vec rscale;
653ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
654ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
655caa4e7f2SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
656caa4e7f2SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
657caa4e7f2SJed Brown       Inject = PETSC_NULL;      /* Callback should create it if it needs Injection */
658caa4e7f2SJed Brown       ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
659caa4e7f2SJed Brown     }
660caa4e7f2SJed Brown   }
661f3fbd535SBarry Smith 
662f3fbd535SBarry Smith   if (!pc->setupcalled) {
663f3fbd535SBarry Smith     for (i=0; i<n; i++) {
664f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
665f3fbd535SBarry Smith     }
666f3fbd535SBarry Smith     for (i=1; i<n; i++) {
667f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
668f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
669f3fbd535SBarry Smith       }
670f3fbd535SBarry Smith     }
671f3fbd535SBarry Smith     for (i=1; i<n; i++) {
672c88c5224SJed Brown       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
673c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
674f3fbd535SBarry Smith     }
675f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
676f3fbd535SBarry Smith       if (!mglevels[i]->b) {
677f3fbd535SBarry Smith         Vec *vec;
678f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
679f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
6806bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
681f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
682f3fbd535SBarry Smith       }
683f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
684f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
685f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
6866bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
687f3fbd535SBarry Smith       }
688f3fbd535SBarry Smith       if (!mglevels[i]->x) {
689f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
690f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
6916bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
692f3fbd535SBarry Smith       }
693f3fbd535SBarry Smith     }
694f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
695f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
696f3fbd535SBarry Smith       Vec *vec;
697f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
698f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
6996bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
700f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
701f3fbd535SBarry Smith     }
702f3fbd535SBarry Smith   }
703f3fbd535SBarry Smith 
70498e04984SBarry Smith   if (pc->dm) {
70598e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
70698e04984SBarry Smith     for (i=0; i<n-1; i++){
70798e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
70898e04984SBarry Smith     }
70998e04984SBarry Smith   }
710f3fbd535SBarry Smith 
711f3fbd535SBarry Smith   for (i=1; i<n; i++) {
712f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
713f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
714f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
715f3fbd535SBarry Smith     }
71663e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
717f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
71863e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
719d42688cbSBarry Smith     if (!mglevels[i]->residual) {
720d42688cbSBarry Smith       Mat mat;
721d42688cbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
722d42688cbSBarry Smith       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
723d42688cbSBarry Smith     }
724f3fbd535SBarry Smith   }
725f3fbd535SBarry Smith   for (i=1; i<n; i++) {
726f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
727f3fbd535SBarry Smith       Mat          downmat,downpmat;
728f3fbd535SBarry Smith       MatStructure matflag;
729f3fbd535SBarry Smith 
730f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
731f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
732f3fbd535SBarry Smith       if (!opsset) {
733f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
734f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
735f3fbd535SBarry Smith       }
736f3fbd535SBarry Smith 
737f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
73863e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
739f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
74063e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
741f3fbd535SBarry Smith     }
742f3fbd535SBarry Smith   }
743f3fbd535SBarry Smith 
744f3fbd535SBarry Smith   /*
745f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
746f3fbd535SBarry Smith   */
747251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
748f3fbd535SBarry Smith   if (preonly) {
749251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
750251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
751251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
752251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
753fd1303e9SJungho Lee     if (!lu && !redundant && !cholesky && !svd) {
754f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
755f3fbd535SBarry Smith     }
756f3fbd535SBarry Smith   }
757f3fbd535SBarry Smith 
758f3fbd535SBarry Smith   if (!pc->setupcalled) {
759f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
760f3fbd535SBarry Smith   }
761f3fbd535SBarry Smith 
76263e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
763f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
76463e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
765f3fbd535SBarry Smith 
766f3fbd535SBarry Smith   /*
767f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
768e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
769f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
770f3fbd535SBarry Smith 
771f3fbd535SBarry Smith    Only support one or the other at the same time.
772f3fbd535SBarry Smith   */
773f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
774acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
775f3fbd535SBarry Smith   if (dump) {
776f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
777f3fbd535SBarry Smith   }
778f3fbd535SBarry Smith   dump = PETSC_FALSE;
779f3fbd535SBarry Smith #endif
780acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
781f3fbd535SBarry Smith   if (dump) {
78232c0f0efSBarry Smith 
783f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
784f3fbd535SBarry Smith   }
785f3fbd535SBarry Smith 
786f3fbd535SBarry Smith   if (viewer) {
787f3fbd535SBarry Smith     for (i=1; i<n; i++) {
788f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
789f3fbd535SBarry Smith     }
790f3fbd535SBarry Smith     for (i=0; i<n; i++) {
791f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
792f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
793f3fbd535SBarry Smith     }
794f3fbd535SBarry Smith   }
795f3fbd535SBarry Smith   PetscFunctionReturn(0);
796f3fbd535SBarry Smith }
797f3fbd535SBarry Smith 
798f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
799f3fbd535SBarry Smith 
800f3fbd535SBarry Smith #undef __FUNCT__
8019dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
8024b9ad928SBarry Smith /*@
80397177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
8044b9ad928SBarry Smith 
8054b9ad928SBarry Smith    Not Collective
8064b9ad928SBarry Smith 
8074b9ad928SBarry Smith    Input Parameter:
8084b9ad928SBarry Smith .  pc - the preconditioner context
8094b9ad928SBarry Smith 
8104b9ad928SBarry Smith    Output parameter:
8114b9ad928SBarry Smith .  levels - the number of levels
8124b9ad928SBarry Smith 
8134b9ad928SBarry Smith    Level: advanced
8144b9ad928SBarry Smith 
8154b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8164b9ad928SBarry Smith 
81797177400SBarry Smith .seealso: PCMGSetLevels()
8184b9ad928SBarry Smith @*/
8197087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8204b9ad928SBarry Smith {
821f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8224b9ad928SBarry Smith 
8234b9ad928SBarry Smith   PetscFunctionBegin;
8240700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8254482741eSBarry Smith   PetscValidIntPointer(levels,2);
826f3fbd535SBarry Smith   *levels = mg->nlevels;
8274b9ad928SBarry Smith   PetscFunctionReturn(0);
8284b9ad928SBarry Smith }
8294b9ad928SBarry Smith 
8304b9ad928SBarry Smith #undef __FUNCT__
8319dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
8324b9ad928SBarry Smith /*@
83397177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8344b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8354b9ad928SBarry Smith 
836ad4df100SBarry Smith    Logically Collective on PC
8374b9ad928SBarry Smith 
8384b9ad928SBarry Smith    Input Parameters:
8394b9ad928SBarry Smith +  pc - the preconditioner context
8409dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
8419dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
8424b9ad928SBarry Smith 
8434b9ad928SBarry Smith    Options Database Key:
8444b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
8454b9ad928SBarry Smith    additive, full, kaskade
8464b9ad928SBarry Smith 
8474b9ad928SBarry Smith    Level: advanced
8484b9ad928SBarry Smith 
8494b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
8504b9ad928SBarry Smith 
85197177400SBarry Smith .seealso: PCMGSetLevels()
8524b9ad928SBarry Smith @*/
8537087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
8544b9ad928SBarry Smith {
855f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
8564b9ad928SBarry Smith 
8574b9ad928SBarry Smith   PetscFunctionBegin;
8580700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
859c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
86031567311SBarry Smith   mg->am = form;
8619dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
8624b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
8634b9ad928SBarry Smith   PetscFunctionReturn(0);
8644b9ad928SBarry Smith }
8654b9ad928SBarry Smith 
8664b9ad928SBarry Smith #undef __FUNCT__
8670d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
8684b9ad928SBarry Smith /*@
8690d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
8704b9ad928SBarry Smith    complicated cycling.
8714b9ad928SBarry Smith 
872ad4df100SBarry Smith    Logically Collective on PC
8734b9ad928SBarry Smith 
8744b9ad928SBarry Smith    Input Parameters:
875c2be2410SBarry Smith +  pc - the multigrid context
8760d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
8774b9ad928SBarry Smith 
8784b9ad928SBarry Smith    Options Database Key:
8790d353602SBarry Smith $  -pc_mg_cycle_type v or w
8804b9ad928SBarry Smith 
8814b9ad928SBarry Smith    Level: advanced
8824b9ad928SBarry Smith 
8834b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8844b9ad928SBarry Smith 
8850d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
8864b9ad928SBarry Smith @*/
8877087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
8884b9ad928SBarry Smith {
889f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
890f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
89179416396SBarry Smith   PetscInt     i,levels;
8924b9ad928SBarry Smith 
8934b9ad928SBarry Smith   PetscFunctionBegin;
8940700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
895e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
896c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
897f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8984b9ad928SBarry Smith 
8994b9ad928SBarry Smith   for (i=0; i<levels; i++) {
900f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
9014b9ad928SBarry Smith   }
9024b9ad928SBarry Smith   PetscFunctionReturn(0);
9034b9ad928SBarry Smith }
9044b9ad928SBarry Smith 
9054b9ad928SBarry Smith #undef __FUNCT__
9068cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
9078cc2d5dfSBarry Smith /*@
9088cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
9098cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9108cc2d5dfSBarry Smith 
911ad4df100SBarry Smith    Logically Collective on PC
9128cc2d5dfSBarry Smith 
9138cc2d5dfSBarry Smith    Input Parameters:
9148cc2d5dfSBarry Smith +  pc - the multigrid context
9158cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
9168cc2d5dfSBarry Smith 
9178cc2d5dfSBarry Smith    Options Database Key:
9188cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
9198cc2d5dfSBarry Smith 
9208cc2d5dfSBarry Smith    Level: advanced
9218cc2d5dfSBarry Smith 
9228cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
9238cc2d5dfSBarry Smith 
9248cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9258cc2d5dfSBarry Smith 
9268cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
9278cc2d5dfSBarry Smith @*/
9287087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
9298cc2d5dfSBarry Smith {
930f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
931f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
9328cc2d5dfSBarry Smith   PetscInt     i,levels;
9338cc2d5dfSBarry Smith 
9348cc2d5dfSBarry Smith   PetscFunctionBegin;
9350700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
936e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
937c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
938f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9398cc2d5dfSBarry Smith 
9408cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
94131567311SBarry Smith     mg->cyclesperpcapply  = n;
9428cc2d5dfSBarry Smith   }
9438cc2d5dfSBarry Smith   PetscFunctionReturn(0);
9448cc2d5dfSBarry Smith }
9458cc2d5dfSBarry Smith 
9468cc2d5dfSBarry Smith #undef __FUNCT__
9479dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
948c2be2410SBarry Smith /*@
94997177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
950c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
951c2be2410SBarry Smith 
952ad4df100SBarry Smith    Logically Collective on PC
953c2be2410SBarry Smith 
954c2be2410SBarry Smith    Input Parameters:
955c91913d3SJed Brown +  pc - the multigrid context
956c91913d3SJed Brown -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
957c2be2410SBarry Smith 
958c2be2410SBarry Smith    Options Database Key:
959c2be2410SBarry Smith $  -pc_mg_galerkin
960c2be2410SBarry Smith 
961c2be2410SBarry Smith    Level: intermediate
962c2be2410SBarry Smith 
963c2be2410SBarry Smith .keywords: MG, set, Galerkin
964c2be2410SBarry Smith 
9653fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
9663fc8bf9cSBarry Smith 
967c2be2410SBarry Smith @*/
968c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
969c2be2410SBarry Smith {
970f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
971c2be2410SBarry Smith 
972c2be2410SBarry Smith   PetscFunctionBegin;
9730700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
974789726d7SBarry Smith   mg->galerkin = use ? 1 : 0;
975c2be2410SBarry Smith   PetscFunctionReturn(0);
976c2be2410SBarry Smith }
977c2be2410SBarry Smith 
978c2be2410SBarry Smith #undef __FUNCT__
9793fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
9803fc8bf9cSBarry Smith /*@
9813fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
9823fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
9833fc8bf9cSBarry Smith 
9843fc8bf9cSBarry Smith    Not Collective
9853fc8bf9cSBarry Smith 
9863fc8bf9cSBarry Smith    Input Parameter:
9873fc8bf9cSBarry Smith .  pc - the multigrid context
9883fc8bf9cSBarry Smith 
9893fc8bf9cSBarry Smith    Output Parameter:
9903fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
9913fc8bf9cSBarry Smith 
9923fc8bf9cSBarry Smith    Options Database Key:
9933fc8bf9cSBarry Smith $  -pc_mg_galerkin
9943fc8bf9cSBarry Smith 
9953fc8bf9cSBarry Smith    Level: intermediate
9963fc8bf9cSBarry Smith 
9973fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
9983fc8bf9cSBarry Smith 
9993fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
10003fc8bf9cSBarry Smith 
10013fc8bf9cSBarry Smith @*/
10027087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
10033fc8bf9cSBarry Smith {
1004f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
10053fc8bf9cSBarry Smith 
10063fc8bf9cSBarry Smith   PetscFunctionBegin;
10070700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
100813fdb3acSJose Roman   *galerkin = (PetscBool)mg->galerkin;
10093fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10103fc8bf9cSBarry Smith }
10113fc8bf9cSBarry Smith 
10123fc8bf9cSBarry Smith #undef __FUNCT__
10139dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
10144b9ad928SBarry Smith /*@
101597177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
101697177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
10174b9ad928SBarry Smith    pre-smoothing steps on different levels.
10184b9ad928SBarry Smith 
1019ad4df100SBarry Smith    Logically Collective on PC
10204b9ad928SBarry Smith 
10214b9ad928SBarry Smith    Input Parameters:
10224b9ad928SBarry Smith +  mg - the multigrid context
10234b9ad928SBarry Smith -  n - the number of smoothing steps
10244b9ad928SBarry Smith 
10254b9ad928SBarry Smith    Options Database Key:
10264b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
10274b9ad928SBarry Smith 
10284b9ad928SBarry Smith    Level: advanced
10294b9ad928SBarry Smith 
10304b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
10314b9ad928SBarry Smith 
103297177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
10334b9ad928SBarry Smith @*/
10347087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
10354b9ad928SBarry Smith {
1036f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1037f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10386849ba73SBarry Smith   PetscErrorCode ierr;
103979416396SBarry Smith   PetscInt       i,levels;
10404b9ad928SBarry Smith 
10414b9ad928SBarry Smith   PetscFunctionBegin;
10420700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1043e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1044c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1045f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10464b9ad928SBarry Smith 
1047b05257ddSBarry Smith   for (i=1; i<levels; i++) {
10484b9ad928SBarry Smith     /* make sure smoother up and down are different */
104997177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1050f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
105131567311SBarry Smith     mg->default_smoothd = n;
10524b9ad928SBarry Smith   }
10534b9ad928SBarry Smith   PetscFunctionReturn(0);
10544b9ad928SBarry Smith }
10554b9ad928SBarry Smith 
10564b9ad928SBarry Smith #undef __FUNCT__
10579dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
10584b9ad928SBarry Smith /*@
105997177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
106097177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
10614b9ad928SBarry Smith    post-smoothing steps on different levels.
10624b9ad928SBarry Smith 
1063ad4df100SBarry Smith    Logically Collective on PC
10644b9ad928SBarry Smith 
10654b9ad928SBarry Smith    Input Parameters:
10664b9ad928SBarry Smith +  mg - the multigrid context
10674b9ad928SBarry Smith -  n - the number of smoothing steps
10684b9ad928SBarry Smith 
10694b9ad928SBarry Smith    Options Database Key:
10704b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
10714b9ad928SBarry Smith 
10724b9ad928SBarry Smith    Level: advanced
10734b9ad928SBarry Smith 
10744b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
1075a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
10764b9ad928SBarry Smith 
10774b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
10784b9ad928SBarry Smith 
107997177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
10804b9ad928SBarry Smith @*/
10817087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
10824b9ad928SBarry Smith {
1083f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1084f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10856849ba73SBarry Smith   PetscErrorCode ierr;
108679416396SBarry Smith   PetscInt       i,levels;
10874b9ad928SBarry Smith 
10884b9ad928SBarry Smith   PetscFunctionBegin;
10890700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1090e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1091c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1092f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10934b9ad928SBarry Smith 
10944b9ad928SBarry Smith   for (i=1; i<levels; i++) {
10954b9ad928SBarry Smith     /* make sure smoother up and down are different */
109697177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1097f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
109831567311SBarry Smith     mg->default_smoothu = n;
10994b9ad928SBarry Smith   }
11004b9ad928SBarry Smith   PetscFunctionReturn(0);
11014b9ad928SBarry Smith }
11024b9ad928SBarry Smith 
11034b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
11044b9ad928SBarry Smith 
11053b09bd56SBarry Smith /*MC
1106ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
11073b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
11083b09bd56SBarry Smith 
11093b09bd56SBarry Smith    Options Database Keys:
11103b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
11110d353602SBarry Smith .  -pc_mg_cycles v or w
111279416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
11133b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
11148c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
11153b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
11163b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
11178cf6037eSBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
11188cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
11198cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1120e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
11218cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
11228cf6037eSBarry Smith                         to the binary output file called binaryoutput
11233b09bd56SBarry Smith 
112424c3aa18SBarry 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
11253b09bd56SBarry Smith 
11268cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
11278cf6037eSBarry Smith 
11283b09bd56SBarry Smith    Level: intermediate
11293b09bd56SBarry Smith 
11308f87f92bSBarry Smith    Concepts: multigrid/multilevel
11313b09bd56SBarry Smith 
11328cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
11330d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
113497177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
113597177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
11360d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11373b09bd56SBarry Smith M*/
11383b09bd56SBarry Smith 
11394b9ad928SBarry Smith EXTERN_C_BEGIN
11404b9ad928SBarry Smith #undef __FUNCT__
11414b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
11427087cfbeSBarry Smith PetscErrorCode  PCCreate_MG(PC pc)
11434b9ad928SBarry Smith {
1144f3fbd535SBarry Smith   PC_MG          *mg;
1145f3fbd535SBarry Smith   PetscErrorCode ierr;
1146f3fbd535SBarry Smith 
11474b9ad928SBarry Smith   PetscFunctionBegin;
1148f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1149f3fbd535SBarry Smith   pc->data    = (void*)mg;
1150f3fbd535SBarry Smith   mg->nlevels = -1;
1151f3fbd535SBarry Smith 
11524b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
11534b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1154a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
11554b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
11564b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
11574b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
11584b9ad928SBarry Smith   PetscFunctionReturn(0);
11594b9ad928SBarry Smith }
11604b9ad928SBarry Smith EXTERN_C_END
1161