xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision d1a3e677c5d75c3517e8677de6882e4d10b79e78)
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;
4525b0b0462SBarry Smith   PetscBool      iascii,isbinary;
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);
457f3fbd535SBarry Smith   if (iascii) {
45831567311SBarry 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);
45931567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
46031567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
461f3fbd535SBarry Smith     }
462218a07d4SBarry Smith     if (mg->galerkin) {
463f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4644f66f45eSBarry Smith     } else {
4654f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
466f3fbd535SBarry Smith     }
467f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
468f3fbd535SBarry Smith       if (!i) {
46989cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
470f3fbd535SBarry Smith       } else {
47189cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
472f3fbd535SBarry Smith       }
473f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
474f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
475f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
476f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
477f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
478f3fbd535SBarry Smith       } else if (i){
47989cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
480f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
481f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
482f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
483f3fbd535SBarry Smith       }
484f3fbd535SBarry Smith     }
4855b0b0462SBarry Smith   } else if (isbinary) {
4865b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
4875b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
4885b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
4895b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
4905b0b0462SBarry Smith       }
4915b0b0462SBarry Smith     }
4925b0b0462SBarry Smith   }
493f3fbd535SBarry Smith   PetscFunctionReturn(0);
494f3fbd535SBarry Smith }
495f3fbd535SBarry Smith 
496b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>
497b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
498cab2e9ccSBarry Smith 
499f3fbd535SBarry Smith /*
500f3fbd535SBarry Smith     Calls setup for the KSP on each level
501f3fbd535SBarry Smith */
502f3fbd535SBarry Smith #undef __FUNCT__
503f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
504f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
505f3fbd535SBarry Smith {
506f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
507f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
508f3fbd535SBarry Smith   PetscErrorCode          ierr;
509f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
51098e04984SBarry Smith   PC                      cpc;
511649052a6SBarry Smith   PetscBool               preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset;
512f3fbd535SBarry Smith   Mat                     dA,dB;
513f3fbd535SBarry Smith   MatStructure            uflag;
514f3fbd535SBarry Smith   Vec                     tvec;
515218a07d4SBarry Smith   DM                      *dms;
516649052a6SBarry Smith   PetscViewer             viewer = 0;
517f3fbd535SBarry Smith 
518f3fbd535SBarry Smith   PetscFunctionBegin;
51901bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5203aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5213aeaf226SBarry Smith     PetscInt levels;
5223aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5233aeaf226SBarry Smith     levels++;
5243aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5253aeaf226SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
5263aeaf226SBarry Smith       n    = levels;
5273aeaf226SBarry Smith       ierr = PCSetFromOptions(pc);CHKERRQ(ierr);  /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5283aeaf226SBarry Smith       mglevels =  mg->levels;
5293aeaf226SBarry Smith     }
5303aeaf226SBarry Smith   }
53198e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5323aeaf226SBarry Smith 
533f3fbd535SBarry Smith 
534f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
535f3fbd535SBarry Smith   /* so use those from global PC */
536f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
537f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
53898e04984SBarry Smith   if (opsset) {
53998e04984SBarry Smith     Mat mmat;
54098e04984SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,PETSC_NULL,&mmat,PETSC_NULL);CHKERRQ(ierr);
54198e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
54298e04984SBarry Smith   }
54398e04984SBarry Smith   if (!opsset) {
544f3fbd535SBarry 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);
545f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
546f3fbd535SBarry Smith   }
547f3fbd535SBarry Smith 
54801bc414fSDmitry 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? */
549ce4cda84SJed Brown   if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
5502d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
5512e499ae9SBarry Smith     Mat p;
55273250ac0SBarry Smith     Vec rscale;
553218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
554218a07d4SBarry Smith     dms[n-1] = pc->dm;
555218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
556fe86f630SJed Brown       KSPDM kdm;
5572ee06e3bSJed Brown       ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);
5583c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
5593c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
560fe86f630SJed Brown       ierr = DMKSPGetContextWrite(dms[i],&kdm);CHKERRQ(ierr);
561*d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
562*d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
563fe86f630SJed Brown       kdm->computerhs = PETSC_NULL;
564fe86f630SJed Brown       kdm->rhsctx = PETSC_NULL;
56511629dbeSBarry Smith       ierr = DMSetFunction(dms[i],0);
56624c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
567e727c939SJed Brown 	ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
5682d2b81a6SBarry Smith 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
56900ac8be1SBarry Smith 	if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
57073250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
5716bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
572218a07d4SBarry Smith       }
57324c3aa18SBarry Smith     }
5742d2b81a6SBarry Smith 
575218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
5766bf464f9SBarry Smith       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
577218a07d4SBarry Smith     }
5782d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
579ce4cda84SJed Brown   }
580cab2e9ccSBarry Smith 
581ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
582ce4cda84SJed Brown     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
583cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
584cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
585218a07d4SBarry Smith   }
586218a07d4SBarry Smith 
587c91913d3SJed Brown   if (mg->galerkin == 1) {
588f3fbd535SBarry Smith     Mat B;
589f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
590f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
591f3fbd535SBarry Smith     if (!pc->setupcalled) {
592f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
593f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
594f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
595f3fbd535SBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
596f3fbd535SBarry Smith         dB   = B;
597f3fbd535SBarry Smith       }
598cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
599f3fbd535SBarry Smith     } else {
600f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
601f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
602f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
603f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
604f3fbd535SBarry Smith         dB   = B;
605f3fbd535SBarry Smith       }
606f3fbd535SBarry Smith     }
607ce4cda84SJed Brown   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
608cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
609cab2e9ccSBarry Smith     for (i=n-2;i>-1; i--) {
610c88c5224SJed Brown       Mat R;
611c88c5224SJed Brown       Vec rscale;
612cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
613cab2e9ccSBarry Smith         Vec *vecs;
614cab2e9ccSBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr);
615cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
616cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
617cab2e9ccSBarry Smith       }
618c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
619c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
620c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
621c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
622cab2e9ccSBarry Smith     }
623f3fbd535SBarry Smith   }
624ccceb50cSJed Brown   if (!mg->galerkin && pc->dm) {
625caa4e7f2SJed Brown     for (i=n-2;i>=0; i--) {
626ccceb50cSJed Brown       DM dmfine,dmcoarse;
627caa4e7f2SJed Brown       Mat Restrict,Inject;
628caa4e7f2SJed Brown       Vec rscale;
629ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
630ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
631caa4e7f2SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
632caa4e7f2SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
633caa4e7f2SJed Brown       Inject = PETSC_NULL;      /* Callback should create it if it needs Injection */
634caa4e7f2SJed Brown       ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
635caa4e7f2SJed Brown     }
636caa4e7f2SJed Brown   }
637f3fbd535SBarry Smith 
638f3fbd535SBarry Smith   if (!pc->setupcalled) {
639f3fbd535SBarry Smith     for (i=0; i<n; i++) {
640f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
641f3fbd535SBarry Smith     }
642f3fbd535SBarry Smith     for (i=1; i<n; i++) {
643f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
644f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
645f3fbd535SBarry Smith       }
646f3fbd535SBarry Smith     }
647f3fbd535SBarry Smith     for (i=1; i<n; i++) {
648c88c5224SJed Brown       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
649c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
650f3fbd535SBarry Smith     }
651f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
652f3fbd535SBarry Smith       if (!mglevels[i]->b) {
653f3fbd535SBarry Smith         Vec *vec;
654f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
655f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
6566bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
657f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
658f3fbd535SBarry Smith       }
659f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
660f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
661f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
6626bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
663f3fbd535SBarry Smith       }
664f3fbd535SBarry Smith       if (!mglevels[i]->x) {
665f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
666f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
6676bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
668f3fbd535SBarry Smith       }
669f3fbd535SBarry Smith     }
670f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
671f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
672f3fbd535SBarry Smith       Vec *vec;
673f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
674f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
6756bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
676f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
677f3fbd535SBarry Smith     }
678f3fbd535SBarry Smith   }
679f3fbd535SBarry Smith 
68098e04984SBarry Smith   if (pc->dm) {
68198e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
68298e04984SBarry Smith     for (i=0; i<n-1; i++){
68398e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
68498e04984SBarry Smith     }
68598e04984SBarry Smith   }
686f3fbd535SBarry Smith 
687f3fbd535SBarry Smith   for (i=1; i<n; i++) {
688f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
689f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
690f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
691f3fbd535SBarry Smith     }
69263e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
693f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
69463e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
695d42688cbSBarry Smith     if (!mglevels[i]->residual) {
696d42688cbSBarry Smith       Mat mat;
697d42688cbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
698d42688cbSBarry Smith       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
699d42688cbSBarry Smith     }
700f3fbd535SBarry Smith   }
701f3fbd535SBarry Smith   for (i=1; i<n; i++) {
702f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
703f3fbd535SBarry Smith       Mat          downmat,downpmat;
704f3fbd535SBarry Smith       MatStructure matflag;
705f3fbd535SBarry Smith 
706f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
707f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
708f3fbd535SBarry Smith       if (!opsset) {
709f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
710f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
711f3fbd535SBarry Smith       }
712f3fbd535SBarry Smith 
713f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
71463e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
715f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
71663e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
717f3fbd535SBarry Smith     }
718f3fbd535SBarry Smith   }
719f3fbd535SBarry Smith 
720f3fbd535SBarry Smith   /*
721f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
722f3fbd535SBarry Smith   */
723251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
724f3fbd535SBarry Smith   if (preonly) {
725251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
726251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
727251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
728251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
729fd1303e9SJungho Lee     if (!lu && !redundant && !cholesky && !svd) {
730f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
731f3fbd535SBarry Smith     }
732f3fbd535SBarry Smith   }
733f3fbd535SBarry Smith 
734f3fbd535SBarry Smith   if (!pc->setupcalled) {
735f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
736f3fbd535SBarry Smith   }
737f3fbd535SBarry Smith 
73863e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
739f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
74063e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
741f3fbd535SBarry Smith 
742f3fbd535SBarry Smith   /*
743f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
744e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
745f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
746f3fbd535SBarry Smith 
747f3fbd535SBarry Smith    Only support one or the other at the same time.
748f3fbd535SBarry Smith   */
749f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
750acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
751f3fbd535SBarry Smith   if (dump) {
752f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
753f3fbd535SBarry Smith   }
754f3fbd535SBarry Smith   dump = PETSC_FALSE;
755f3fbd535SBarry Smith #endif
756acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
757f3fbd535SBarry Smith   if (dump) {
75832c0f0efSBarry Smith 
759f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
760f3fbd535SBarry Smith   }
761f3fbd535SBarry Smith 
762f3fbd535SBarry Smith   if (viewer) {
763f3fbd535SBarry Smith     for (i=1; i<n; i++) {
764f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
765f3fbd535SBarry Smith     }
766f3fbd535SBarry Smith     for (i=0; i<n; i++) {
767f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
768f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
769f3fbd535SBarry Smith     }
770f3fbd535SBarry Smith   }
771f3fbd535SBarry Smith   PetscFunctionReturn(0);
772f3fbd535SBarry Smith }
773f3fbd535SBarry Smith 
774f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
775f3fbd535SBarry Smith 
776f3fbd535SBarry Smith #undef __FUNCT__
7779dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
7784b9ad928SBarry Smith /*@
77997177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
7804b9ad928SBarry Smith 
7814b9ad928SBarry Smith    Not Collective
7824b9ad928SBarry Smith 
7834b9ad928SBarry Smith    Input Parameter:
7844b9ad928SBarry Smith .  pc - the preconditioner context
7854b9ad928SBarry Smith 
7864b9ad928SBarry Smith    Output parameter:
7874b9ad928SBarry Smith .  levels - the number of levels
7884b9ad928SBarry Smith 
7894b9ad928SBarry Smith    Level: advanced
7904b9ad928SBarry Smith 
7914b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
7924b9ad928SBarry Smith 
79397177400SBarry Smith .seealso: PCMGSetLevels()
7944b9ad928SBarry Smith @*/
7957087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
7964b9ad928SBarry Smith {
797f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
7984b9ad928SBarry Smith 
7994b9ad928SBarry Smith   PetscFunctionBegin;
8000700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8014482741eSBarry Smith   PetscValidIntPointer(levels,2);
802f3fbd535SBarry Smith   *levels = mg->nlevels;
8034b9ad928SBarry Smith   PetscFunctionReturn(0);
8044b9ad928SBarry Smith }
8054b9ad928SBarry Smith 
8064b9ad928SBarry Smith #undef __FUNCT__
8079dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
8084b9ad928SBarry Smith /*@
80997177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8104b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8114b9ad928SBarry Smith 
812ad4df100SBarry Smith    Logically Collective on PC
8134b9ad928SBarry Smith 
8144b9ad928SBarry Smith    Input Parameters:
8154b9ad928SBarry Smith +  pc - the preconditioner context
8169dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
8179dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
8184b9ad928SBarry Smith 
8194b9ad928SBarry Smith    Options Database Key:
8204b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
8214b9ad928SBarry Smith    additive, full, kaskade
8224b9ad928SBarry Smith 
8234b9ad928SBarry Smith    Level: advanced
8244b9ad928SBarry Smith 
8254b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
8264b9ad928SBarry Smith 
82797177400SBarry Smith .seealso: PCMGSetLevels()
8284b9ad928SBarry Smith @*/
8297087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
8304b9ad928SBarry Smith {
831f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
8324b9ad928SBarry Smith 
8334b9ad928SBarry Smith   PetscFunctionBegin;
8340700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
835c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
83631567311SBarry Smith   mg->am = form;
8379dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
8384b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
8394b9ad928SBarry Smith   PetscFunctionReturn(0);
8404b9ad928SBarry Smith }
8414b9ad928SBarry Smith 
8424b9ad928SBarry Smith #undef __FUNCT__
8430d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
8444b9ad928SBarry Smith /*@
8450d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
8464b9ad928SBarry Smith    complicated cycling.
8474b9ad928SBarry Smith 
848ad4df100SBarry Smith    Logically Collective on PC
8494b9ad928SBarry Smith 
8504b9ad928SBarry Smith    Input Parameters:
851c2be2410SBarry Smith +  pc - the multigrid context
8520d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
8534b9ad928SBarry Smith 
8544b9ad928SBarry Smith    Options Database Key:
8550d353602SBarry Smith $  -pc_mg_cycle_type v or w
8564b9ad928SBarry Smith 
8574b9ad928SBarry Smith    Level: advanced
8584b9ad928SBarry Smith 
8594b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8604b9ad928SBarry Smith 
8610d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
8624b9ad928SBarry Smith @*/
8637087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
8644b9ad928SBarry Smith {
865f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
866f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
86779416396SBarry Smith   PetscInt     i,levels;
8684b9ad928SBarry Smith 
8694b9ad928SBarry Smith   PetscFunctionBegin;
8700700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
871e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
872c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
873f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8744b9ad928SBarry Smith 
8754b9ad928SBarry Smith   for (i=0; i<levels; i++) {
876f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
8774b9ad928SBarry Smith   }
8784b9ad928SBarry Smith   PetscFunctionReturn(0);
8794b9ad928SBarry Smith }
8804b9ad928SBarry Smith 
8814b9ad928SBarry Smith #undef __FUNCT__
8828cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
8838cc2d5dfSBarry Smith /*@
8848cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
8858cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
8868cc2d5dfSBarry Smith 
887ad4df100SBarry Smith    Logically Collective on PC
8888cc2d5dfSBarry Smith 
8898cc2d5dfSBarry Smith    Input Parameters:
8908cc2d5dfSBarry Smith +  pc - the multigrid context
8918cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
8928cc2d5dfSBarry Smith 
8938cc2d5dfSBarry Smith    Options Database Key:
8948cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
8958cc2d5dfSBarry Smith 
8968cc2d5dfSBarry Smith    Level: advanced
8978cc2d5dfSBarry Smith 
8988cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
8998cc2d5dfSBarry Smith 
9008cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9018cc2d5dfSBarry Smith 
9028cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
9038cc2d5dfSBarry Smith @*/
9047087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
9058cc2d5dfSBarry Smith {
906f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
907f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
9088cc2d5dfSBarry Smith   PetscInt     i,levels;
9098cc2d5dfSBarry Smith 
9108cc2d5dfSBarry Smith   PetscFunctionBegin;
9110700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
912e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
913c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
914f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9158cc2d5dfSBarry Smith 
9168cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
91731567311SBarry Smith     mg->cyclesperpcapply  = n;
9188cc2d5dfSBarry Smith   }
9198cc2d5dfSBarry Smith   PetscFunctionReturn(0);
9208cc2d5dfSBarry Smith }
9218cc2d5dfSBarry Smith 
9228cc2d5dfSBarry Smith #undef __FUNCT__
9239dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
924c2be2410SBarry Smith /*@
92597177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
926c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
927c2be2410SBarry Smith 
928ad4df100SBarry Smith    Logically Collective on PC
929c2be2410SBarry Smith 
930c2be2410SBarry Smith    Input Parameters:
931c91913d3SJed Brown +  pc - the multigrid context
932c91913d3SJed Brown -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
933c2be2410SBarry Smith 
934c2be2410SBarry Smith    Options Database Key:
935c2be2410SBarry Smith $  -pc_mg_galerkin
936c2be2410SBarry Smith 
937c2be2410SBarry Smith    Level: intermediate
938c2be2410SBarry Smith 
939c2be2410SBarry Smith .keywords: MG, set, Galerkin
940c2be2410SBarry Smith 
9413fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
9423fc8bf9cSBarry Smith 
943c2be2410SBarry Smith @*/
944c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
945c2be2410SBarry Smith {
946f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
947c2be2410SBarry Smith 
948c2be2410SBarry Smith   PetscFunctionBegin;
9490700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
950789726d7SBarry Smith   mg->galerkin = use ? 1 : 0;
951c2be2410SBarry Smith   PetscFunctionReturn(0);
952c2be2410SBarry Smith }
953c2be2410SBarry Smith 
954c2be2410SBarry Smith #undef __FUNCT__
9553fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
9563fc8bf9cSBarry Smith /*@
9573fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
9583fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
9593fc8bf9cSBarry Smith 
9603fc8bf9cSBarry Smith    Not Collective
9613fc8bf9cSBarry Smith 
9623fc8bf9cSBarry Smith    Input Parameter:
9633fc8bf9cSBarry Smith .  pc - the multigrid context
9643fc8bf9cSBarry Smith 
9653fc8bf9cSBarry Smith    Output Parameter:
9663fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
9673fc8bf9cSBarry Smith 
9683fc8bf9cSBarry Smith    Options Database Key:
9693fc8bf9cSBarry Smith $  -pc_mg_galerkin
9703fc8bf9cSBarry Smith 
9713fc8bf9cSBarry Smith    Level: intermediate
9723fc8bf9cSBarry Smith 
9733fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
9743fc8bf9cSBarry Smith 
9753fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
9763fc8bf9cSBarry Smith 
9773fc8bf9cSBarry Smith @*/
9787087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
9793fc8bf9cSBarry Smith {
980f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
9813fc8bf9cSBarry Smith 
9823fc8bf9cSBarry Smith   PetscFunctionBegin;
9830700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
98413fdb3acSJose Roman   *galerkin = (PetscBool)mg->galerkin;
9853fc8bf9cSBarry Smith   PetscFunctionReturn(0);
9863fc8bf9cSBarry Smith }
9873fc8bf9cSBarry Smith 
9883fc8bf9cSBarry Smith #undef __FUNCT__
9899dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
9904b9ad928SBarry Smith /*@
99197177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
99297177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
9934b9ad928SBarry Smith    pre-smoothing steps on different levels.
9944b9ad928SBarry Smith 
995ad4df100SBarry Smith    Logically Collective on PC
9964b9ad928SBarry Smith 
9974b9ad928SBarry Smith    Input Parameters:
9984b9ad928SBarry Smith +  mg - the multigrid context
9994b9ad928SBarry Smith -  n - the number of smoothing steps
10004b9ad928SBarry Smith 
10014b9ad928SBarry Smith    Options Database Key:
10024b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
10034b9ad928SBarry Smith 
10044b9ad928SBarry Smith    Level: advanced
10054b9ad928SBarry Smith 
10064b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
10074b9ad928SBarry Smith 
100897177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
10094b9ad928SBarry Smith @*/
10107087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
10114b9ad928SBarry Smith {
1012f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1013f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10146849ba73SBarry Smith   PetscErrorCode ierr;
101579416396SBarry Smith   PetscInt       i,levels;
10164b9ad928SBarry Smith 
10174b9ad928SBarry Smith   PetscFunctionBegin;
10180700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1019e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1020c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1021f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10224b9ad928SBarry Smith 
1023b05257ddSBarry Smith   for (i=1; i<levels; i++) {
10244b9ad928SBarry Smith     /* make sure smoother up and down are different */
102597177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1026f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
102731567311SBarry Smith     mg->default_smoothd = n;
10284b9ad928SBarry Smith   }
10294b9ad928SBarry Smith   PetscFunctionReturn(0);
10304b9ad928SBarry Smith }
10314b9ad928SBarry Smith 
10324b9ad928SBarry Smith #undef __FUNCT__
10339dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
10344b9ad928SBarry Smith /*@
103597177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
103697177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
10374b9ad928SBarry Smith    post-smoothing steps on different levels.
10384b9ad928SBarry Smith 
1039ad4df100SBarry Smith    Logically Collective on PC
10404b9ad928SBarry Smith 
10414b9ad928SBarry Smith    Input Parameters:
10424b9ad928SBarry Smith +  mg - the multigrid context
10434b9ad928SBarry Smith -  n - the number of smoothing steps
10444b9ad928SBarry Smith 
10454b9ad928SBarry Smith    Options Database Key:
10464b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
10474b9ad928SBarry Smith 
10484b9ad928SBarry Smith    Level: advanced
10494b9ad928SBarry Smith 
10504b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
1051a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
10524b9ad928SBarry Smith 
10534b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
10544b9ad928SBarry Smith 
105597177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
10564b9ad928SBarry Smith @*/
10577087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
10584b9ad928SBarry Smith {
1059f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1060f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10616849ba73SBarry Smith   PetscErrorCode ierr;
106279416396SBarry Smith   PetscInt       i,levels;
10634b9ad928SBarry Smith 
10644b9ad928SBarry Smith   PetscFunctionBegin;
10650700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1066e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1067c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1068f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10694b9ad928SBarry Smith 
10704b9ad928SBarry Smith   for (i=1; i<levels; i++) {
10714b9ad928SBarry Smith     /* make sure smoother up and down are different */
107297177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1073f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
107431567311SBarry Smith     mg->default_smoothu = n;
10754b9ad928SBarry Smith   }
10764b9ad928SBarry Smith   PetscFunctionReturn(0);
10774b9ad928SBarry Smith }
10784b9ad928SBarry Smith 
10794b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
10804b9ad928SBarry Smith 
10813b09bd56SBarry Smith /*MC
1082ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
10833b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
10843b09bd56SBarry Smith 
10853b09bd56SBarry Smith    Options Database Keys:
10863b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
10870d353602SBarry Smith .  -pc_mg_cycles v or w
108879416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
10893b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
10908c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
10913b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
10923b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
10938cf6037eSBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
10948cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
10958cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1096e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
10978cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
10988cf6037eSBarry Smith                         to the binary output file called binaryoutput
10993b09bd56SBarry Smith 
110024c3aa18SBarry 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
11013b09bd56SBarry Smith 
11028cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
11038cf6037eSBarry Smith 
11043b09bd56SBarry Smith    Level: intermediate
11053b09bd56SBarry Smith 
11068f87f92bSBarry Smith    Concepts: multigrid/multilevel
11073b09bd56SBarry Smith 
11088cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
11090d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
111097177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
111197177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
11120d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11133b09bd56SBarry Smith M*/
11143b09bd56SBarry Smith 
11154b9ad928SBarry Smith EXTERN_C_BEGIN
11164b9ad928SBarry Smith #undef __FUNCT__
11174b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
11187087cfbeSBarry Smith PetscErrorCode  PCCreate_MG(PC pc)
11194b9ad928SBarry Smith {
1120f3fbd535SBarry Smith   PC_MG          *mg;
1121f3fbd535SBarry Smith   PetscErrorCode ierr;
1122f3fbd535SBarry Smith 
11234b9ad928SBarry Smith   PetscFunctionBegin;
1124f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1125f3fbd535SBarry Smith   pc->data    = (void*)mg;
1126f3fbd535SBarry Smith   mg->nlevels = -1;
1127f3fbd535SBarry Smith 
11284b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
11294b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1130a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
11314b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
11324b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
11334b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
11344b9ad928SBarry Smith   PetscFunctionReturn(0);
11354b9ad928SBarry Smith }
11364b9ad928SBarry Smith EXTERN_C_END
1137