xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 3aeaf22671e2e6825835bd6d8cd7fec9c5bd6fd9)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscpcmg.h" I*/
64b9ad928SBarry Smith 
74b9ad928SBarry Smith 
84b9ad928SBarry Smith #undef __FUNCT__
99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private"
1031567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
114b9ad928SBarry Smith {
1231567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1331567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
146849ba73SBarry Smith   PetscErrorCode ierr;
1531567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
184b9ad928SBarry Smith 
1963e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2031567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
2163e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2231567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2363e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2431567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2563e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
264b9ad928SBarry Smith 
274b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
2831567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
294b9ad928SBarry Smith       PetscReal rnorm;
3031567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
314b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3270441072SBarry Smith         if (rnorm < mg->abstol) {
334d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
341e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
354b9ad928SBarry Smith         } else {
364d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
371e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr);
384b9ad928SBarry Smith         }
394b9ad928SBarry Smith         PetscFunctionReturn(0);
404b9ad928SBarry Smith       }
414b9ad928SBarry Smith     }
424b9ad928SBarry Smith 
4331567311SBarry Smith     mgc = *(mglevelsin - 1);
4463e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4531567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
4663e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
47efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
484b9ad928SBarry Smith     while (cycles--) {
4931567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
504b9ad928SBarry Smith     }
5163e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5231567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5363e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5463e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5531567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
5663e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith   }
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
614b9ad928SBarry Smith #undef __FUNCT__
624b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
63ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool  zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
644b9ad928SBarry Smith {
65f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
66f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
67dfbe8321SBarry Smith   PetscErrorCode ierr;
68f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
694b9ad928SBarry Smith 
704b9ad928SBarry Smith   PetscFunctionBegin;
71f3fbd535SBarry Smith   mglevels[levels-1]->b    = b;
72f3fbd535SBarry Smith   mglevels[levels-1]->x    = x;
734b9ad928SBarry Smith 
7431567311SBarry Smith   mg->rtol = rtol;
7531567311SBarry Smith   mg->abstol = abstol;
7631567311SBarry Smith   mg->dtol = dtol;
774b9ad928SBarry Smith   if (rtol) {
784b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
794b9ad928SBarry Smith     PetscReal rnorm;
807319c654SBarry Smith     if (zeroguess) {
817319c654SBarry Smith       ierr               = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
827319c654SBarry Smith     } else {
83f3fbd535SBarry Smith       ierr               = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
844b9ad928SBarry Smith       ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
857319c654SBarry Smith     }
8631567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
8770441072SBarry Smith   } else if (abstol) {
8831567311SBarry Smith     mg->ttol = abstol;
894b9ad928SBarry Smith   } else {
9031567311SBarry Smith     mg->ttol = 0.0;
914b9ad928SBarry Smith   }
924b9ad928SBarry Smith 
934d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
944d0a8057SBarry Smith      stop prematurely do to small residual */
954d0a8057SBarry Smith   for (i=1; i<levels; i++) {
96f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
97f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
98f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
994b9ad928SBarry Smith     }
1004d0a8057SBarry Smith   }
1014d0a8057SBarry Smith 
1024d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1034d0a8057SBarry Smith   for (i=0; i<its; i++) {
104f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1054d0a8057SBarry Smith     if (*reason) break;
1064d0a8057SBarry Smith   }
1074d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1084d0a8057SBarry Smith   *outits = i;
1094b9ad928SBarry Smith   PetscFunctionReturn(0);
1104b9ad928SBarry Smith }
1114b9ad928SBarry Smith 
1124b9ad928SBarry Smith #undef __FUNCT__
113*3aeaf226SBarry Smith #define __FUNCT__ "PCReset_MG"
114*3aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
115*3aeaf226SBarry Smith {
116*3aeaf226SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
117*3aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
118*3aeaf226SBarry Smith   PetscErrorCode ierr;
119*3aeaf226SBarry Smith   PetscInt       i,n;
120*3aeaf226SBarry Smith 
121*3aeaf226SBarry Smith   PetscFunctionBegin;
122*3aeaf226SBarry Smith   if (mglevels) {
123*3aeaf226SBarry Smith     n = mglevels[0]->levels;
124*3aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
125*3aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
126*3aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
127*3aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
128*3aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
129*3aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
130*3aeaf226SBarry Smith     }
131*3aeaf226SBarry Smith 
132*3aeaf226SBarry Smith     for (i=0; i<n; i++) {
133*3aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
134*3aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
135*3aeaf226SBarry Smith 	ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
136*3aeaf226SBarry Smith       }
137*3aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
138*3aeaf226SBarry Smith     }
139*3aeaf226SBarry Smith   }
140*3aeaf226SBarry Smith   PetscFunctionReturn(0);
141*3aeaf226SBarry Smith }
142*3aeaf226SBarry Smith 
143*3aeaf226SBarry Smith #undef __FUNCT__
1449dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
1454b9ad928SBarry Smith /*@C
14697177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1474b9ad928SBarry Smith    Must be called before any other MG routine.
1484b9ad928SBarry Smith 
149ad4df100SBarry Smith    Logically Collective on PC
1504b9ad928SBarry Smith 
1514b9ad928SBarry Smith    Input Parameters:
1524b9ad928SBarry Smith +  pc - the preconditioner context
1534b9ad928SBarry Smith .  levels - the number of levels
1544b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1554b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith    Level: intermediate
1584b9ad928SBarry Smith 
1594b9ad928SBarry Smith    Notes:
1604b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1614b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1624b9ad928SBarry Smith 
1634b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1644b9ad928SBarry Smith 
16597177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1664b9ad928SBarry Smith @*/
1677087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1684b9ad928SBarry Smith {
169dfbe8321SBarry Smith   PetscErrorCode ierr;
170f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
171f3fbd535SBarry Smith   MPI_Comm       comm = ((PetscObject)pc)->comm;
172*3aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
173f3fbd535SBarry Smith   PetscInt       i;
174f3fbd535SBarry Smith   PetscMPIInt    size;
175f3fbd535SBarry Smith   const char     *prefix;
176f3fbd535SBarry Smith   PC             ipc;
177*3aeaf226SBarry Smith   PetscInt       n;
1784b9ad928SBarry Smith 
1794b9ad928SBarry Smith   PetscFunctionBegin;
1800700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
181548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
182548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
183*3aeaf226SBarry Smith   if (mglevels) {
184*3aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
185*3aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
186*3aeaf226SBarry Smith     n = mglevels[0]->levels;
187*3aeaf226SBarry Smith     for (i=0; i<n; i++) {
188*3aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
189*3aeaf226SBarry Smith 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
190*3aeaf226SBarry Smith       }
191*3aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
192*3aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
193*3aeaf226SBarry Smith     }
194*3aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
195*3aeaf226SBarry Smith   }
196f3fbd535SBarry Smith 
197f3fbd535SBarry Smith   mg->nlevels      = levels;
198218a07d4SBarry Smith   mg->galerkin     = PETSC_FALSE;
199218a07d4SBarry Smith   mg->galerkinused = PETSC_FALSE;
200f3fbd535SBarry Smith 
201f3fbd535SBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr);
202f3fbd535SBarry Smith   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
203f3fbd535SBarry Smith 
204f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
205f3fbd535SBarry Smith 
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;
21131567311SBarry Smith     mg->default_smoothu = 1;
21231567311SBarry Smith     mg->default_smoothd = 1;
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);
220f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
22131567311SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
222f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
223f3fbd535SBarry Smith 
224f3fbd535SBarry Smith     /* do special stuff for coarse grid */
225f3fbd535SBarry Smith     if (!i && levels > 1) {
226f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
227f3fbd535SBarry Smith 
228f3fbd535SBarry Smith       /* coarse solve is (redundant) LU by default */
229f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
230f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
231f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
232f3fbd535SBarry Smith       if (size > 1) {
233f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
234f3fbd535SBarry Smith       } else {
235f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
236f3fbd535SBarry Smith       }
237f3fbd535SBarry Smith 
238f3fbd535SBarry Smith     } else {
239f3fbd535SBarry Smith       char tprefix[128];
240f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
241f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
242f3fbd535SBarry Smith     }
243f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr);
244f3fbd535SBarry Smith     mglevels[i]->smoothu    = mglevels[i]->smoothd;
24531567311SBarry Smith     mg->rtol                = 0.0;
24631567311SBarry Smith     mg->abstol              = 0.0;
24731567311SBarry Smith     mg->dtol                = 0.0;
24831567311SBarry Smith     mg->ttol                = 0.0;
24931567311SBarry Smith     mg->cyclesperpcapply    = 1;
250f3fbd535SBarry Smith   }
25131567311SBarry Smith   mg->am          = PC_MG_MULTIPLICATIVE;
252f3fbd535SBarry Smith   mg->levels      = mglevels;
2534b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
2544b9ad928SBarry Smith   PetscFunctionReturn(0);
2554b9ad928SBarry Smith }
2564b9ad928SBarry Smith 
257c07bf074SBarry Smith 
258c07bf074SBarry Smith #undef __FUNCT__
259c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
260c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
261c07bf074SBarry Smith {
262c07bf074SBarry Smith   PetscErrorCode ierr;
263a06653b4SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
264a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
265a06653b4SBarry Smith   PetscInt       i,n;
266c07bf074SBarry Smith 
267c07bf074SBarry Smith   PetscFunctionBegin;
268a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
269a06653b4SBarry Smith   if (mglevels) {
270a06653b4SBarry Smith     n = mglevels[0]->levels;
271a06653b4SBarry Smith     for (i=0; i<n; i++) {
272a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2736bf464f9SBarry Smith 	ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
274a06653b4SBarry Smith       }
2756bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
276a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
277a06653b4SBarry Smith     }
278a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
279a06653b4SBarry Smith   }
280c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
281f3fbd535SBarry Smith   PetscFunctionReturn(0);
282f3fbd535SBarry Smith }
283f3fbd535SBarry Smith 
284f3fbd535SBarry Smith 
285f3fbd535SBarry Smith 
28609573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
28709573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
28809573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
289f3fbd535SBarry Smith 
290f3fbd535SBarry Smith /*
291f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
292f3fbd535SBarry Smith              or full cycle of multigrid.
293f3fbd535SBarry Smith 
294f3fbd535SBarry Smith   Note:
295f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
296f3fbd535SBarry Smith */
297f3fbd535SBarry Smith #undef __FUNCT__
298f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
299f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
300f3fbd535SBarry Smith {
301f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
302f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
303f3fbd535SBarry Smith   PetscErrorCode ierr;
304f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
305f3fbd535SBarry Smith 
306f3fbd535SBarry Smith   PetscFunctionBegin;
307e1d8e5deSBarry Smith 
308e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
309298cc208SBarry Smith   for (i=0; i<levels; i++) {
310e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
311e1d8e5deSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
312298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
313e1d8e5deSBarry Smith     }
314e1d8e5deSBarry Smith   }
315e1d8e5deSBarry Smith 
316f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
317f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
31831567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
319f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
32031567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
321f3fbd535SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
322f3fbd535SBarry Smith     }
323f3fbd535SBarry Smith   }
32431567311SBarry Smith   else if (mg->am == PC_MG_ADDITIVE) {
32531567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
326f3fbd535SBarry Smith   }
32731567311SBarry Smith   else if (mg->am == PC_MG_KASKADE) {
32831567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
329f3fbd535SBarry Smith   }
330f3fbd535SBarry Smith   else {
331f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
332f3fbd535SBarry Smith   }
333f3fbd535SBarry Smith   PetscFunctionReturn(0);
334f3fbd535SBarry Smith }
335f3fbd535SBarry Smith 
336f3fbd535SBarry Smith 
337f3fbd535SBarry Smith #undef __FUNCT__
338f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
339f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc)
340f3fbd535SBarry Smith {
341f3fbd535SBarry Smith   PetscErrorCode ierr;
342f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
343ace3abfcSBarry Smith   PetscBool      flg;
344f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
345f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
346f3fbd535SBarry Smith   PCMGType       mgtype;
347f3fbd535SBarry Smith   PCMGCycleType  mgctype;
348f3fbd535SBarry Smith 
349f3fbd535SBarry Smith   PetscFunctionBegin;
350f3fbd535SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
351*3aeaf226SBarry Smith     if (!mg->levels) {
352f3fbd535SBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
353eb3f98d2SBarry Smith       if (!flg && pc->dm) {
354eb3f98d2SBarry Smith         ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
355eb3f98d2SBarry Smith         levels++;
356*3aeaf226SBarry Smith         mg->usedmfornumberoflevels = PETSC_TRUE;
357eb3f98d2SBarry Smith       }
358f3fbd535SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
359f3fbd535SBarry Smith     }
360*3aeaf226SBarry Smith     mglevels = mg->levels;
361*3aeaf226SBarry Smith 
362f3fbd535SBarry Smith     mgctype = (PCMGCycleType) mglevels[0]->cycles;
363f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
364f3fbd535SBarry Smith     if (flg) {
365f3fbd535SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
366f3fbd535SBarry Smith     };
367f3fbd535SBarry Smith     flg  = PETSC_FALSE;
368acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
369f3fbd535SBarry Smith     if (flg) {
370f3fbd535SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
371f3fbd535SBarry Smith     }
372f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
373f3fbd535SBarry Smith     if (flg) {
374f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
375f3fbd535SBarry Smith     }
376f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
377f3fbd535SBarry Smith     if (flg) {
378f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
379f3fbd535SBarry Smith     }
38031567311SBarry Smith     mgtype = mg->am;
381f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
382f3fbd535SBarry Smith     if (flg) {
383f3fbd535SBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
384f3fbd535SBarry Smith     }
38531567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
38631567311SBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
387f3fbd535SBarry Smith       if (flg) {
388f3fbd535SBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
389f3fbd535SBarry Smith       }
390f3fbd535SBarry Smith     }
391f3fbd535SBarry Smith     flg  = PETSC_FALSE;
392acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
393f3fbd535SBarry Smith     if (flg) {
394f3fbd535SBarry Smith       PetscInt i;
395f3fbd535SBarry Smith       char     eventname[128];
39663e6d426SJed Brown       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
397f3fbd535SBarry Smith       levels = mglevels[0]->levels;
398f3fbd535SBarry Smith       for (i=0; i<levels; i++) {
399f3fbd535SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
40063e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
401f3fbd535SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
40263e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
403f3fbd535SBarry Smith         if (i) {
404f3fbd535SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
40563e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
406f3fbd535SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
40763e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
408f3fbd535SBarry Smith         }
409f3fbd535SBarry Smith       }
410f3fbd535SBarry Smith     }
411f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
412f3fbd535SBarry Smith   PetscFunctionReturn(0);
413f3fbd535SBarry Smith }
414f3fbd535SBarry Smith 
415f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
416f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
417f3fbd535SBarry Smith 
418f3fbd535SBarry Smith #undef __FUNCT__
419f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
420f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
421f3fbd535SBarry Smith {
422f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
423f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
424f3fbd535SBarry Smith   PetscErrorCode ierr;
425f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
426ace3abfcSBarry Smith   PetscBool      iascii;
427f3fbd535SBarry Smith 
428f3fbd535SBarry Smith   PetscFunctionBegin;
4292692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
430f3fbd535SBarry Smith   if (iascii) {
43131567311SBarry 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);
43231567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
43331567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
434f3fbd535SBarry Smith     }
435218a07d4SBarry Smith     if (mg->galerkin) {
436f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4374f66f45eSBarry Smith     } else {
4384f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
439f3fbd535SBarry Smith     }
440f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
441f3fbd535SBarry Smith       if (!i) {
4422d19a89fSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
443f3fbd535SBarry Smith       } else {
44431567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr);
445f3fbd535SBarry Smith       }
446f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
447f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
448f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
449f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
450f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
451f3fbd535SBarry Smith       } else if (i){
45231567311SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr);
453f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
454f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
455f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
456f3fbd535SBarry Smith       }
457f3fbd535SBarry Smith     }
458f3fbd535SBarry Smith   } else {
45965e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
460f3fbd535SBarry Smith   }
461f3fbd535SBarry Smith   PetscFunctionReturn(0);
462f3fbd535SBarry Smith }
463f3fbd535SBarry Smith 
464f3fbd535SBarry Smith /*
465f3fbd535SBarry Smith     Calls setup for the KSP on each level
466f3fbd535SBarry Smith */
467f3fbd535SBarry Smith #undef __FUNCT__
468f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
469f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
470f3fbd535SBarry Smith {
471f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
472f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
473f3fbd535SBarry Smith   PetscErrorCode          ierr;
474f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
475f3fbd535SBarry Smith   PC                      cpc,mpc;
476ace3abfcSBarry Smith   PetscBool               preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset;
477f3fbd535SBarry Smith   PetscViewerASCIIMonitor ascii;
478f3fbd535SBarry Smith   PetscViewer             viewer = PETSC_NULL;
479f3fbd535SBarry Smith   MPI_Comm                comm;
480f3fbd535SBarry Smith   Mat                     dA,dB;
481f3fbd535SBarry Smith   MatStructure            uflag;
482f3fbd535SBarry Smith   Vec                     tvec;
483218a07d4SBarry Smith   DM                      *dms;
484f3fbd535SBarry Smith 
485f3fbd535SBarry Smith   PetscFunctionBegin;
486*3aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
487*3aeaf226SBarry Smith     PetscInt levels;
488*3aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
489*3aeaf226SBarry Smith     levels++;
490*3aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
491*3aeaf226SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
492*3aeaf226SBarry Smith       n    = levels;
493*3aeaf226SBarry Smith       ierr = PCSetFromOptions(pc);CHKERRQ(ierr);  /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
494*3aeaf226SBarry Smith       mglevels =  mg->levels;
495*3aeaf226SBarry Smith     }
496*3aeaf226SBarry Smith   }
497*3aeaf226SBarry Smith 
498f3fbd535SBarry Smith 
499f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
500f3fbd535SBarry Smith   /* so use those from global PC */
501f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
502f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
503f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
504f3fbd535SBarry Smith   ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr);
5051338a6b9SJed Brown   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) {
506f3fbd535SBarry 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);
507f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
508f3fbd535SBarry Smith   }
509f3fbd535SBarry Smith 
5102d2b81a6SBarry Smith   if (pc->dm && !pc->setupcalled) {
5112d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
5122e499ae9SBarry Smith     Mat p;
513218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
514218a07d4SBarry Smith     dms[n-1] = pc->dm;
515218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
516218a07d4SBarry Smith       ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr);
5173c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
5183c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
51911629dbeSBarry Smith       ierr = DMSetFunction(dms[i],0);
52011629dbeSBarry Smith       ierr = DMSetInitialGuess(dms[i],0);
52124c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
5222d2b81a6SBarry Smith 	ierr = DMGetInterpolation(dms[i],dms[i+1],&p,PETSC_NULL);CHKERRQ(ierr);
5232d2b81a6SBarry Smith 	ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
5246bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
525218a07d4SBarry Smith       }
52624c3aa18SBarry Smith     }
5272d2b81a6SBarry Smith 
528218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
5296bf464f9SBarry Smith       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
530218a07d4SBarry Smith     }
5312d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
532218a07d4SBarry Smith   }
533218a07d4SBarry Smith 
534218a07d4SBarry Smith   if (mg->galerkin) {
535f3fbd535SBarry Smith     Mat B;
536218a07d4SBarry Smith     mg->galerkinused = PETSC_TRUE;
537f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
538f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
539f3fbd535SBarry Smith     if (!pc->setupcalled) {
540f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
541f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
542f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
543f3fbd535SBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
544f3fbd535SBarry Smith         dB   = B;
545f3fbd535SBarry Smith       }
546cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
547f3fbd535SBarry Smith     } else {
548f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
549f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
550f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
551f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
552f3fbd535SBarry Smith         dB   = B;
553f3fbd535SBarry Smith       }
554f3fbd535SBarry Smith     }
555f3fbd535SBarry Smith   }
556f3fbd535SBarry Smith 
557f3fbd535SBarry Smith   if (!pc->setupcalled) {
558acfcf0e5SJed Brown     ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr);
559f3fbd535SBarry Smith 
560f3fbd535SBarry Smith     for (i=0; i<n; i++) {
561f3fbd535SBarry Smith       if (monitor) {
562f3fbd535SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr);
563f3fbd535SBarry Smith         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
564c2efdce3SBarry Smith         ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
565f3fbd535SBarry Smith       }
566f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
567f3fbd535SBarry Smith     }
568f3fbd535SBarry Smith     for (i=1; i<n; i++) {
569f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
570f3fbd535SBarry Smith         if (monitor) {
571f3fbd535SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr);
572f3fbd535SBarry Smith           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
573c2efdce3SBarry Smith           ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
574f3fbd535SBarry Smith         }
575f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
576f3fbd535SBarry Smith       }
577f3fbd535SBarry Smith     }
578f3fbd535SBarry Smith     for (i=1; i<n; i++) {
579f3fbd535SBarry Smith       if (mglevels[i]->restrct && !mglevels[i]->interpolate) {
580f3fbd535SBarry Smith         ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr);
581f3fbd535SBarry Smith       }
582f3fbd535SBarry Smith       if (!mglevels[i]->restrct && mglevels[i]->interpolate) {
583f3fbd535SBarry Smith         ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr);
584f3fbd535SBarry Smith       }
585f3fbd535SBarry Smith #if defined(PETSC_USE_DEBUG)
586f3fbd535SBarry Smith       if (!mglevels[i]->restrct || !mglevels[i]->interpolate) {
58765e19b50SBarry Smith         SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
588f3fbd535SBarry Smith       }
589f3fbd535SBarry Smith #endif
590f3fbd535SBarry Smith     }
591f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
592f3fbd535SBarry Smith       if (!mglevels[i]->b) {
593f3fbd535SBarry Smith         Vec *vec;
594f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
595f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
5966bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
597f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
598f3fbd535SBarry Smith       }
599f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
600f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
601f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
6026bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
603f3fbd535SBarry Smith       }
604f3fbd535SBarry Smith       if (!mglevels[i]->x) {
605f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
606f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
6076bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
608f3fbd535SBarry Smith       }
609f3fbd535SBarry Smith     }
610f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
611f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
612f3fbd535SBarry Smith       Vec *vec;
613f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
614f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
6156bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
616f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
617f3fbd535SBarry Smith     }
618f3fbd535SBarry Smith   }
619f3fbd535SBarry Smith 
620f3fbd535SBarry Smith 
621f3fbd535SBarry Smith   for (i=1; i<n; i++) {
622f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
623f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
624f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
625f3fbd535SBarry Smith     }
62663e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
627f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
62863e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
629d42688cbSBarry Smith     if (!mglevels[i]->residual) {
630d42688cbSBarry Smith       Mat mat;
631d42688cbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
632d42688cbSBarry Smith       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
633d42688cbSBarry Smith     }
634f3fbd535SBarry Smith   }
635f3fbd535SBarry Smith   for (i=1; i<n; i++) {
636f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
637f3fbd535SBarry Smith       Mat          downmat,downpmat;
638f3fbd535SBarry Smith       MatStructure matflag;
639ace3abfcSBarry Smith       PetscBool    opsset;
640f3fbd535SBarry Smith 
641f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
642f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
643f3fbd535SBarry Smith       if (!opsset) {
644f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
645f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
646f3fbd535SBarry Smith       }
647f3fbd535SBarry Smith 
648f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
64963e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
650f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
65163e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
652f3fbd535SBarry Smith     }
653f3fbd535SBarry Smith   }
654f3fbd535SBarry Smith 
655f3fbd535SBarry Smith   /*
656f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
657f3fbd535SBarry Smith   */
658f3fbd535SBarry Smith   ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
659f3fbd535SBarry Smith   if (preonly) {
660f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
661f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
662f3fbd535SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
663f3fbd535SBarry Smith     if (!lu && !redundant && !cholesky) {
664f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
665f3fbd535SBarry Smith     }
666f3fbd535SBarry Smith   }
667f3fbd535SBarry Smith 
668f3fbd535SBarry Smith   if (!pc->setupcalled) {
669f3fbd535SBarry Smith     if (monitor) {
670f3fbd535SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr);
671f3fbd535SBarry Smith       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
672c2efdce3SBarry Smith       ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void**))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
673f3fbd535SBarry Smith     }
674f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
675f3fbd535SBarry Smith   }
676f3fbd535SBarry Smith 
67763e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
678f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
67963e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
680f3fbd535SBarry Smith 
681f3fbd535SBarry Smith   /*
682f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
683e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
684f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
685f3fbd535SBarry Smith 
686f3fbd535SBarry Smith    Only support one or the other at the same time.
687f3fbd535SBarry Smith   */
688f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
689acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
690f3fbd535SBarry Smith   if (dump) {
691f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
692f3fbd535SBarry Smith   }
693f3fbd535SBarry Smith   dump = PETSC_FALSE;
694f3fbd535SBarry Smith #endif
695acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
696f3fbd535SBarry Smith   if (dump) {
697f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
698f3fbd535SBarry Smith   }
699f3fbd535SBarry Smith 
700f3fbd535SBarry Smith   if (viewer) {
701f3fbd535SBarry Smith     for (i=1; i<n; i++) {
702f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
703f3fbd535SBarry Smith     }
704f3fbd535SBarry Smith     for (i=0; i<n; i++) {
705f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
706f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
707f3fbd535SBarry Smith     }
708f3fbd535SBarry Smith   }
709f3fbd535SBarry Smith   PetscFunctionReturn(0);
710f3fbd535SBarry Smith }
711f3fbd535SBarry Smith 
712f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
713f3fbd535SBarry Smith 
714f3fbd535SBarry Smith #undef __FUNCT__
7159dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
7164b9ad928SBarry Smith /*@
71797177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
7184b9ad928SBarry Smith 
7194b9ad928SBarry Smith    Not Collective
7204b9ad928SBarry Smith 
7214b9ad928SBarry Smith    Input Parameter:
7224b9ad928SBarry Smith .  pc - the preconditioner context
7234b9ad928SBarry Smith 
7244b9ad928SBarry Smith    Output parameter:
7254b9ad928SBarry Smith .  levels - the number of levels
7264b9ad928SBarry Smith 
7274b9ad928SBarry Smith    Level: advanced
7284b9ad928SBarry Smith 
7294b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
7304b9ad928SBarry Smith 
73197177400SBarry Smith .seealso: PCMGSetLevels()
7324b9ad928SBarry Smith @*/
7337087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
7344b9ad928SBarry Smith {
735f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
7364b9ad928SBarry Smith 
7374b9ad928SBarry Smith   PetscFunctionBegin;
7380700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
7394482741eSBarry Smith   PetscValidIntPointer(levels,2);
740f3fbd535SBarry Smith   *levels = mg->nlevels;
7414b9ad928SBarry Smith   PetscFunctionReturn(0);
7424b9ad928SBarry Smith }
7434b9ad928SBarry Smith 
7444b9ad928SBarry Smith #undef __FUNCT__
7459dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
7464b9ad928SBarry Smith /*@
74797177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
7484b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
7494b9ad928SBarry Smith 
750ad4df100SBarry Smith    Logically Collective on PC
7514b9ad928SBarry Smith 
7524b9ad928SBarry Smith    Input Parameters:
7534b9ad928SBarry Smith +  pc - the preconditioner context
7549dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
7559dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
7564b9ad928SBarry Smith 
7574b9ad928SBarry Smith    Options Database Key:
7584b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
7594b9ad928SBarry Smith    additive, full, kaskade
7604b9ad928SBarry Smith 
7614b9ad928SBarry Smith    Level: advanced
7624b9ad928SBarry Smith 
7634b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
7644b9ad928SBarry Smith 
76597177400SBarry Smith .seealso: PCMGSetLevels()
7664b9ad928SBarry Smith @*/
7677087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
7684b9ad928SBarry Smith {
769f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
7704b9ad928SBarry Smith 
7714b9ad928SBarry Smith   PetscFunctionBegin;
7720700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
773c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
77431567311SBarry Smith   mg->am = form;
7759dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
7764b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
7774b9ad928SBarry Smith   PetscFunctionReturn(0);
7784b9ad928SBarry Smith }
7794b9ad928SBarry Smith 
7804b9ad928SBarry Smith #undef __FUNCT__
7810d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
7824b9ad928SBarry Smith /*@
7830d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
7844b9ad928SBarry Smith    complicated cycling.
7854b9ad928SBarry Smith 
786ad4df100SBarry Smith    Logically Collective on PC
7874b9ad928SBarry Smith 
7884b9ad928SBarry Smith    Input Parameters:
789c2be2410SBarry Smith +  pc - the multigrid context
7900d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
7914b9ad928SBarry Smith 
7924b9ad928SBarry Smith    Options Database Key:
7930d353602SBarry Smith $  -pc_mg_cycle_type v or w
7944b9ad928SBarry Smith 
7954b9ad928SBarry Smith    Level: advanced
7964b9ad928SBarry Smith 
7974b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7984b9ad928SBarry Smith 
7990d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
8004b9ad928SBarry Smith @*/
8017087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
8024b9ad928SBarry Smith {
803f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
804f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
80579416396SBarry Smith   PetscInt     i,levels;
8064b9ad928SBarry Smith 
8074b9ad928SBarry Smith   PetscFunctionBegin;
8080700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
809e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
810c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
811f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8124b9ad928SBarry Smith 
8134b9ad928SBarry Smith   for (i=0; i<levels; i++) {
814f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
8154b9ad928SBarry Smith   }
8164b9ad928SBarry Smith   PetscFunctionReturn(0);
8174b9ad928SBarry Smith }
8184b9ad928SBarry Smith 
8194b9ad928SBarry Smith #undef __FUNCT__
8208cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
8218cc2d5dfSBarry Smith /*@
8228cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
8238cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
8248cc2d5dfSBarry Smith 
825ad4df100SBarry Smith    Logically Collective on PC
8268cc2d5dfSBarry Smith 
8278cc2d5dfSBarry Smith    Input Parameters:
8288cc2d5dfSBarry Smith +  pc - the multigrid context
8298cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
8308cc2d5dfSBarry Smith 
8318cc2d5dfSBarry Smith    Options Database Key:
8328cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
8338cc2d5dfSBarry Smith 
8348cc2d5dfSBarry Smith    Level: advanced
8358cc2d5dfSBarry Smith 
8368cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
8378cc2d5dfSBarry Smith 
8388cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8398cc2d5dfSBarry Smith 
8408cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
8418cc2d5dfSBarry Smith @*/
8427087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
8438cc2d5dfSBarry Smith {
844f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
845f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
8468cc2d5dfSBarry Smith   PetscInt     i,levels;
8478cc2d5dfSBarry Smith 
8488cc2d5dfSBarry Smith   PetscFunctionBegin;
8490700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
850e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
851c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
852f3fbd535SBarry Smith   levels = mglevels[0]->levels;
8538cc2d5dfSBarry Smith 
8548cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
85531567311SBarry Smith     mg->cyclesperpcapply  = n;
8568cc2d5dfSBarry Smith   }
8578cc2d5dfSBarry Smith   PetscFunctionReturn(0);
8588cc2d5dfSBarry Smith }
8598cc2d5dfSBarry Smith 
8608cc2d5dfSBarry Smith #undef __FUNCT__
8619dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
862c2be2410SBarry Smith /*@
86397177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
864c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
865c2be2410SBarry Smith 
866ad4df100SBarry Smith    Logically Collective on PC
867c2be2410SBarry Smith 
868c2be2410SBarry Smith    Input Parameters:
8693fc8bf9cSBarry Smith .  pc - the multigrid context
870c2be2410SBarry Smith 
871c2be2410SBarry Smith    Options Database Key:
872c2be2410SBarry Smith $  -pc_mg_galerkin
873c2be2410SBarry Smith 
874c2be2410SBarry Smith    Level: intermediate
875c2be2410SBarry Smith 
876c2be2410SBarry Smith .keywords: MG, set, Galerkin
877c2be2410SBarry Smith 
8783fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
8793fc8bf9cSBarry Smith 
880c2be2410SBarry Smith @*/
8817087cfbeSBarry Smith PetscErrorCode  PCMGSetGalerkin(PC pc)
882c2be2410SBarry Smith {
883f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
884c2be2410SBarry Smith 
885c2be2410SBarry Smith   PetscFunctionBegin;
8860700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
887218a07d4SBarry Smith   mg->galerkin = PETSC_TRUE;
888c2be2410SBarry Smith   PetscFunctionReturn(0);
889c2be2410SBarry Smith }
890c2be2410SBarry Smith 
891c2be2410SBarry Smith #undef __FUNCT__
8923fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
8933fc8bf9cSBarry Smith /*@
8943fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
8953fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
8963fc8bf9cSBarry Smith 
8973fc8bf9cSBarry Smith    Not Collective
8983fc8bf9cSBarry Smith 
8993fc8bf9cSBarry Smith    Input Parameter:
9003fc8bf9cSBarry Smith .  pc - the multigrid context
9013fc8bf9cSBarry Smith 
9023fc8bf9cSBarry Smith    Output Parameter:
9033fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
9043fc8bf9cSBarry Smith 
9053fc8bf9cSBarry Smith    Options Database Key:
9063fc8bf9cSBarry Smith $  -pc_mg_galerkin
9073fc8bf9cSBarry Smith 
9083fc8bf9cSBarry Smith    Level: intermediate
9093fc8bf9cSBarry Smith 
9103fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
9113fc8bf9cSBarry Smith 
9123fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
9133fc8bf9cSBarry Smith 
9143fc8bf9cSBarry Smith @*/
9157087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
9163fc8bf9cSBarry Smith {
917f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
9183fc8bf9cSBarry Smith 
9193fc8bf9cSBarry Smith   PetscFunctionBegin;
9200700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
921218a07d4SBarry Smith   *galerkin = mg->galerkin;
9223fc8bf9cSBarry Smith   PetscFunctionReturn(0);
9233fc8bf9cSBarry Smith }
9243fc8bf9cSBarry Smith 
9253fc8bf9cSBarry Smith #undef __FUNCT__
9269dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
9274b9ad928SBarry Smith /*@
92897177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
92997177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
9304b9ad928SBarry Smith    pre-smoothing steps on different levels.
9314b9ad928SBarry Smith 
932ad4df100SBarry Smith    Logically Collective on PC
9334b9ad928SBarry Smith 
9344b9ad928SBarry Smith    Input Parameters:
9354b9ad928SBarry Smith +  mg - the multigrid context
9364b9ad928SBarry Smith -  n - the number of smoothing steps
9374b9ad928SBarry Smith 
9384b9ad928SBarry Smith    Options Database Key:
9394b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
9404b9ad928SBarry Smith 
9414b9ad928SBarry Smith    Level: advanced
9424b9ad928SBarry Smith 
9434b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
9444b9ad928SBarry Smith 
94597177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
9464b9ad928SBarry Smith @*/
9477087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
9484b9ad928SBarry Smith {
949f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
950f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9516849ba73SBarry Smith   PetscErrorCode ierr;
95279416396SBarry Smith   PetscInt       i,levels;
9534b9ad928SBarry Smith 
9544b9ad928SBarry Smith   PetscFunctionBegin;
9550700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
956e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
957c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
958f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9594b9ad928SBarry Smith 
960b05257ddSBarry Smith   for (i=1; i<levels; i++) {
9614b9ad928SBarry Smith     /* make sure smoother up and down are different */
96297177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
963f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
96431567311SBarry Smith     mg->default_smoothd = n;
9654b9ad928SBarry Smith   }
9664b9ad928SBarry Smith   PetscFunctionReturn(0);
9674b9ad928SBarry Smith }
9684b9ad928SBarry Smith 
9694b9ad928SBarry Smith #undef __FUNCT__
9709dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
9714b9ad928SBarry Smith /*@
97297177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
97397177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
9744b9ad928SBarry Smith    post-smoothing steps on different levels.
9754b9ad928SBarry Smith 
976ad4df100SBarry Smith    Logically Collective on PC
9774b9ad928SBarry Smith 
9784b9ad928SBarry Smith    Input Parameters:
9794b9ad928SBarry Smith +  mg - the multigrid context
9804b9ad928SBarry Smith -  n - the number of smoothing steps
9814b9ad928SBarry Smith 
9824b9ad928SBarry Smith    Options Database Key:
9834b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
9844b9ad928SBarry Smith 
9854b9ad928SBarry Smith    Level: advanced
9864b9ad928SBarry Smith 
9874b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
988a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
9894b9ad928SBarry Smith 
9904b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
9914b9ad928SBarry Smith 
99297177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
9934b9ad928SBarry Smith @*/
9947087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
9954b9ad928SBarry Smith {
996f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
997f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
9986849ba73SBarry Smith   PetscErrorCode ierr;
99979416396SBarry Smith   PetscInt       i,levels;
10004b9ad928SBarry Smith 
10014b9ad928SBarry Smith   PetscFunctionBegin;
10020700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1003e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1004c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1005f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10064b9ad928SBarry Smith 
10074b9ad928SBarry Smith   for (i=1; i<levels; i++) {
10084b9ad928SBarry Smith     /* make sure smoother up and down are different */
100997177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1010f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
101131567311SBarry Smith     mg->default_smoothu = n;
10124b9ad928SBarry Smith   }
10134b9ad928SBarry Smith   PetscFunctionReturn(0);
10144b9ad928SBarry Smith }
10154b9ad928SBarry Smith 
10164b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
10174b9ad928SBarry Smith 
10183b09bd56SBarry Smith /*MC
1019ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
10203b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
10213b09bd56SBarry Smith 
10223b09bd56SBarry Smith    Options Database Keys:
10233b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
10240d353602SBarry Smith .  -pc_mg_cycles v or w
102579416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
10263b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
10273b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
10283b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
10293b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
103068eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
10313b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1032e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
10333b09bd56SBarry Smith 
103424c3aa18SBarry 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
10353b09bd56SBarry Smith 
10363b09bd56SBarry Smith    Level: intermediate
10373b09bd56SBarry Smith 
10388f87f92bSBarry Smith    Concepts: multigrid/multilevel
10393b09bd56SBarry Smith 
104024c3aa18SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC
10410d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
104297177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
104397177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
10440d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
10453b09bd56SBarry Smith M*/
10463b09bd56SBarry Smith 
10474b9ad928SBarry Smith EXTERN_C_BEGIN
10484b9ad928SBarry Smith #undef __FUNCT__
10494b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
10507087cfbeSBarry Smith PetscErrorCode  PCCreate_MG(PC pc)
10514b9ad928SBarry Smith {
1052f3fbd535SBarry Smith   PC_MG          *mg;
1053f3fbd535SBarry Smith   PetscErrorCode ierr;
1054f3fbd535SBarry Smith 
10554b9ad928SBarry Smith   PetscFunctionBegin;
1056f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1057f3fbd535SBarry Smith   pc->data    = (void*)mg;
1058f3fbd535SBarry Smith   mg->nlevels = -1;
1059f3fbd535SBarry Smith 
10604b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
10614b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1062a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
10634b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
10644b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
10654b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
10664b9ad928SBarry Smith   PetscFunctionReturn(0);
10674b9ad928SBarry Smith }
10684b9ad928SBarry Smith EXTERN_C_END
1069