xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 4d0a8057dc766494e73f985a8c029a3f3edde3ad)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
54b9ad928SBarry Smith */
64b9ad928SBarry Smith #include "src/ksp/pc/impls/mg/mgimpl.h"                    /*I "petscmg.h" I*/
74b9ad928SBarry Smith 
84b9ad928SBarry Smith 
94b9ad928SBarry Smith #undef __FUNCT__
109dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private"
11*4d0a8057SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG **mglevels,PCRichardsonConvergedReason *reason)
124b9ad928SBarry Smith {
139dcbbd2bSBarry Smith   PC_MG          *mg = *mglevels,*mgc;
146849ba73SBarry Smith   PetscErrorCode ierr;
150d353602SBarry Smith   PetscInt       cycles = (PetscInt) mg->cycles;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
184b9ad928SBarry Smith 
1932cf1786SBarry Smith   if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
200d353602SBarry Smith   ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr);  /* pre-smooth */
2132cf1786SBarry Smith   if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
224b9ad928SBarry Smith   if (mg->level) {  /* not the coarsest grid */
2332cf1786SBarry Smith     if (mg->eventresidual) {ierr = PetscLogEventBegin(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
244b9ad928SBarry Smith     ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr);
2532cf1786SBarry Smith     if (mg->eventresidual) {ierr = PetscLogEventEnd(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);}
264b9ad928SBarry Smith 
274b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
28*4d0a8057SBarry Smith     if (mg->level == mg->levels-1 && mg->ttol && reason) {
294b9ad928SBarry Smith       PetscReal rnorm;
304b9ad928SBarry Smith       ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr);
314b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3270441072SBarry Smith         if (rnorm < mg->abstol) {
33*4d0a8057SBarry 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 {
36*4d0a8057SBarry 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 
434b9ad928SBarry Smith     mgc = *(mglevels - 1);
4432cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
454b9ad928SBarry Smith     ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr);
4632cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
47efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
484b9ad928SBarry Smith     while (cycles--) {
49*4d0a8057SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels-1,reason);CHKERRQ(ierr);
504b9ad928SBarry Smith     }
5132cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
524b9ad928SBarry Smith     ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr);
5332cf1786SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5432cf1786SBarry Smith     if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
550d353602SBarry Smith     ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr);    /* post smooth */
5632cf1786SBarry Smith     if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith   }
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
614b9ad928SBarry Smith /*
629dcbbd2bSBarry Smith        PCMGCreate_Private - Creates a PC_MG structure for use with the
634b9ad928SBarry Smith                multigrid code. Level 0 is the coarsest. (But the
644b9ad928SBarry Smith                finest level is stored first in the array).
654b9ad928SBarry Smith 
664b9ad928SBarry Smith */
674b9ad928SBarry Smith #undef __FUNCT__
689dcbbd2bSBarry Smith #define __FUNCT__ "PCMGCreate_Private"
699dcbbd2bSBarry Smith static PetscErrorCode PCMGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,PC_MG ***result)
704b9ad928SBarry Smith {
719dcbbd2bSBarry Smith   PC_MG          **mg;
726849ba73SBarry Smith   PetscErrorCode ierr;
7379416396SBarry Smith   PetscInt       i;
7479416396SBarry Smith   PetscMPIInt    size;
75f69a0ea3SMatthew Knepley   const char     *prefix;
764b9ad928SBarry Smith   PC             ipc;
774b9ad928SBarry Smith 
784b9ad928SBarry Smith   PetscFunctionBegin;
799dcbbd2bSBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mg);CHKERRQ(ierr);
8038f2d2fdSLisandro Dalcin   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
814b9ad928SBarry Smith 
824b9ad928SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
834b9ad928SBarry Smith 
844b9ad928SBarry Smith   for (i=0; i<levels; i++) {
8538f2d2fdSLisandro Dalcin     ierr = PetscNewLog(pc,PC_MG,&mg[i]);CHKERRQ(ierr);
864b9ad928SBarry Smith     mg[i]->level           = i;
874b9ad928SBarry Smith     mg[i]->levels          = levels;
880d353602SBarry Smith     mg[i]->cycles          = PC_MG_CYCLE_V;
89c2be2410SBarry Smith     mg[i]->galerkin        = PETSC_FALSE;
90c2be2410SBarry Smith     mg[i]->galerkinused    = PETSC_FALSE;
9179416396SBarry Smith     mg[i]->default_smoothu = 1;
9279416396SBarry Smith     mg[i]->default_smoothd = 1;
934b9ad928SBarry Smith 
944b9ad928SBarry Smith     if (comms) comm = comms[i];
954b9ad928SBarry Smith     ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr);
9679416396SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr);
974b9ad928SBarry Smith     ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr);
984b9ad928SBarry Smith 
994b9ad928SBarry Smith     /* do special stuff for coarse grid */
1004b9ad928SBarry Smith     if (!i && levels > 1) {
1014b9ad928SBarry Smith       ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
1024b9ad928SBarry Smith 
1034b9ad928SBarry Smith       /* coarse solve is (redundant) LU by default */
1044b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
1054b9ad928SBarry Smith       ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr);
1064b9ad928SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1074b9ad928SBarry Smith       if (size > 1) {
1084b9ad928SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
1090d7810c8SBarry Smith       } else {
1104b9ad928SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
1110d7810c8SBarry Smith       }
1124b9ad928SBarry Smith 
1134b9ad928SBarry Smith     } else {
1142e70eadcSBarry Smith       char tprefix[128];
11513f74950SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
1162e70eadcSBarry Smith       ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr);
1174b9ad928SBarry Smith     }
11852e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr);
1194b9ad928SBarry Smith     mg[i]->smoothu             = mg[i]->smoothd;
1204b9ad928SBarry Smith     mg[i]->rtol                = 0.0;
12170441072SBarry Smith     mg[i]->abstol              = 0.0;
1224b9ad928SBarry Smith     mg[i]->dtol                = 0.0;
1234b9ad928SBarry Smith     mg[i]->ttol                = 0.0;
12432cf1786SBarry Smith     mg[i]->eventsmoothsetup    = 0;
12532cf1786SBarry Smith     mg[i]->eventsmoothsolve    = 0;
12632cf1786SBarry Smith     mg[i]->eventresidual       = 0;
12732cf1786SBarry Smith     mg[i]->eventinterprestrict = 0;
1288cc2d5dfSBarry Smith     mg[i]->cyclesperpcapply    = 1;
1294b9ad928SBarry Smith   }
1304b9ad928SBarry Smith   *result = mg;
1314b9ad928SBarry Smith   PetscFunctionReturn(0);
1324b9ad928SBarry Smith }
1334b9ad928SBarry Smith 
1344b9ad928SBarry Smith #undef __FUNCT__
1354b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_MG"
1366849ba73SBarry Smith static PetscErrorCode PCDestroy_MG(PC pc)
1374b9ad928SBarry Smith {
1389dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1396849ba73SBarry Smith   PetscErrorCode ierr;
14038f2d2fdSLisandro Dalcin   PetscInt       i,n;
1414b9ad928SBarry Smith 
1424b9ad928SBarry Smith   PetscFunctionBegin;
14338f2d2fdSLisandro Dalcin   if (!mg) PetscFunctionReturn(0);
14438f2d2fdSLisandro Dalcin   n = mg[0]->levels;
145fccaa45eSBarry Smith   for (i=0; i<n-1; i++) {
146630ba2eeSBarry Smith     if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);}
147fccaa45eSBarry Smith     if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);}
148fccaa45eSBarry Smith     if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);}
149fccaa45eSBarry Smith     if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);}
150fccaa45eSBarry Smith     if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);}
151fccaa45eSBarry Smith   }
152fccaa45eSBarry Smith 
1534b9ad928SBarry Smith   for (i=0; i<n; i++) {
1544b9ad928SBarry Smith     if (mg[i]->smoothd != mg[i]->smoothu) {
1554b9ad928SBarry Smith       ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr);
1564b9ad928SBarry Smith     }
1574b9ad928SBarry Smith     ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr);
1584b9ad928SBarry Smith     ierr = PetscFree(mg[i]);CHKERRQ(ierr);
1594b9ad928SBarry Smith   }
1604b9ad928SBarry Smith   ierr = PetscFree(mg);CHKERRQ(ierr);
1614b9ad928SBarry Smith   PetscFunctionReturn(0);
1624b9ad928SBarry Smith }
1634b9ad928SBarry Smith 
1644b9ad928SBarry Smith 
1654b9ad928SBarry Smith 
1669dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**);
1671e2582c4SBarry Smith EXTERN PetscErrorCode PCMGFCycle_Private(PC,PC_MG**);
1689dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**);
1694b9ad928SBarry Smith 
1704b9ad928SBarry Smith /*
1714b9ad928SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
1724b9ad928SBarry Smith              or full cycle of multigrid.
1734b9ad928SBarry Smith 
1744b9ad928SBarry Smith   Note:
1759dcbbd2bSBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
1764b9ad928SBarry Smith */
1774b9ad928SBarry Smith #undef __FUNCT__
1784b9ad928SBarry Smith #define __FUNCT__ "PCApply_MG"
1796849ba73SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
1804b9ad928SBarry Smith {
1819dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1826849ba73SBarry Smith   PetscErrorCode ierr;
1838cc2d5dfSBarry Smith   PetscInt       levels = mg[0]->levels,i;
1844b9ad928SBarry Smith 
1854b9ad928SBarry Smith   PetscFunctionBegin;
1864b9ad928SBarry Smith   mg[levels-1]->b = b;
1874b9ad928SBarry Smith   mg[levels-1]->x = x;
188ef70c39aSMatthew Knepley   if (!mg[levels-1]->r && mg[0]->am != PC_MG_ADDITIVE && levels > 1) {
1890a6bb862SBarry Smith     Vec tvec;
1900a6bb862SBarry Smith     ierr = VecDuplicate(mg[levels-1]->b,&tvec);CHKERRQ(ierr);
19197177400SBarry Smith     ierr = PCMGSetR(pc,levels-1,tvec);CHKERRQ(ierr);
1920a6bb862SBarry Smith     ierr = VecDestroy(tvec);CHKERRQ(ierr);
1930a6bb862SBarry Smith   }
1949dcbbd2bSBarry Smith   if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
195efb30889SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1968cc2d5dfSBarry Smith     for (i=0; i<mg[0]->cyclesperpcapply; i++) {
1971e2582c4SBarry Smith       ierr = PCMGMCycle_Private(pc,mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
1984b9ad928SBarry Smith     }
1998cc2d5dfSBarry Smith   }
2009dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_ADDITIVE) {
2019dcbbd2bSBarry Smith     ierr = PCMGACycle_Private(mg);CHKERRQ(ierr);
2024b9ad928SBarry Smith   }
2039dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_KASKADE) {
2049dcbbd2bSBarry Smith     ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr);
2054b9ad928SBarry Smith   }
2064b9ad928SBarry Smith   else {
2071e2582c4SBarry Smith     ierr = PCMGFCycle_Private(pc,mg);CHKERRQ(ierr);
2084b9ad928SBarry Smith   }
2094b9ad928SBarry Smith   PetscFunctionReturn(0);
2104b9ad928SBarry Smith }
2114b9ad928SBarry Smith 
2124b9ad928SBarry Smith #undef __FUNCT__
2134b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
214*4d0a8057SBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscInt *outits,PCRichardsonConvergedReason *reason)
2154b9ad928SBarry Smith {
2169dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
217dfbe8321SBarry Smith   PetscErrorCode ierr;
218*4d0a8057SBarry Smith   PetscInt       levels = mg[0]->levels,i;
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith   PetscFunctionBegin;
2214b9ad928SBarry Smith   mg[levels-1]->b    = b;
2224b9ad928SBarry Smith   mg[levels-1]->x    = x;
2234b9ad928SBarry Smith 
2244b9ad928SBarry Smith   mg[levels-1]->rtol = rtol;
22570441072SBarry Smith   mg[levels-1]->abstol = abstol;
2264b9ad928SBarry Smith   mg[levels-1]->dtol = dtol;
2274b9ad928SBarry Smith   if (rtol) {
2284b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
2294b9ad928SBarry Smith     PetscReal rnorm;
2304b9ad928SBarry Smith     ierr               = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr);
2314b9ad928SBarry Smith     ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
23270441072SBarry Smith     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
23370441072SBarry Smith   } else if (abstol) {
23470441072SBarry Smith     mg[levels-1]->ttol = abstol;
2354b9ad928SBarry Smith   } else {
2364b9ad928SBarry Smith     mg[levels-1]->ttol = 0.0;
2374b9ad928SBarry Smith   }
2384b9ad928SBarry Smith 
239*4d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
240*4d0a8057SBarry Smith      stop prematurely do to small residual */
241*4d0a8057SBarry Smith   for (i=1; i<levels; i++) {
242*4d0a8057SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
243*4d0a8057SBarry Smith     if (mg[i]->smoothu != mg[i]->smoothd) {
244*4d0a8057SBarry Smith       ierr = KSPSetTolerances(mg[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2454b9ad928SBarry Smith     }
246*4d0a8057SBarry Smith   }
247*4d0a8057SBarry Smith 
248*4d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
249*4d0a8057SBarry Smith   for (i=0; i<its; i++) {
250*4d0a8057SBarry Smith     ierr = PCMGMCycle_Private(pc,mg+levels-1,reason);CHKERRQ(ierr);
251*4d0a8057SBarry Smith     if (*reason) break;
252*4d0a8057SBarry Smith   }
253*4d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
254*4d0a8057SBarry Smith   *outits = i;
2554b9ad928SBarry Smith   PetscFunctionReturn(0);
2564b9ad928SBarry Smith }
2574b9ad928SBarry Smith 
2584b9ad928SBarry Smith #undef __FUNCT__
2594b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
2606ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_MG(PC pc)
2614b9ad928SBarry Smith {
262dfbe8321SBarry Smith   PetscErrorCode ierr;
2638cc2d5dfSBarry Smith   PetscInt       m,levels = 1,cycles;
2644b9ad928SBarry Smith   PetscTruth     flg;
2659dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
26690da95b0SMatthew Knepley   PCMGType       mgtype = PC_MG_ADDITIVE;
2671aa34eeaSBarry Smith   PCMGCycleType  mgctype;
2684b9ad928SBarry Smith 
2694b9ad928SBarry Smith   PetscFunctionBegin;
2704b9ad928SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
2714b9ad928SBarry Smith     if (!pc->data) {
2729dcbbd2bSBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
27397177400SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
274cf502942SBarry Smith       mg = (PC_MG**)pc->data;
2754b9ad928SBarry Smith     }
276f2070a76SMatthew Knepley     mgctype = (PCMGCycleType) mg[0]->cycles;
2770d353602SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
2784b9ad928SBarry Smith     if (flg) {
2790d353602SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
2800d353602SBarry Smith     };
2819dcbbd2bSBarry Smith     ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr);
282c2be2410SBarry Smith     if (flg) {
28397177400SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
284c2be2410SBarry Smith     }
2859dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
2864b9ad928SBarry Smith     if (flg) {
28797177400SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
2884b9ad928SBarry Smith     }
2899dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
2904b9ad928SBarry Smith     if (flg) {
29197177400SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
2924b9ad928SBarry Smith     }
2939dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
2948cc2d5dfSBarry Smith     if (flg) {
2958cc2d5dfSBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2968cc2d5dfSBarry Smith     }
2978cc2d5dfSBarry Smith     if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
2988cc2d5dfSBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
2998cc2d5dfSBarry Smith       if (flg) {
3008cc2d5dfSBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
3018cc2d5dfSBarry Smith       }
3028cc2d5dfSBarry Smith     }
3034b9ad928SBarry Smith     ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr);
3044b9ad928SBarry Smith     if (flg) {
3054f5ab15aSBarry Smith       PetscInt i;
3064b9ad928SBarry Smith       char     eventname[128];
3074b9ad928SBarry Smith       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
3084b9ad928SBarry Smith       levels = mg[0]->levels;
3094b9ad928SBarry Smith       for (i=0; i<levels; i++) {
31032cf1786SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
3118cbcd9ccSBarry Smith         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventsmoothsetup);CHKERRQ(ierr);
31232cf1786SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
3138cbcd9ccSBarry Smith         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventsmoothsolve);CHKERRQ(ierr);
31432cf1786SBarry Smith         if (i) {
31532cf1786SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
316a3bc4eb9SBarry Smith           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventresidual);CHKERRQ(ierr);
31732cf1786SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
3188cbcd9ccSBarry Smith           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventinterprestrict);CHKERRQ(ierr);
31932cf1786SBarry Smith         }
3204b9ad928SBarry Smith       }
3214b9ad928SBarry Smith     }
3224b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3234b9ad928SBarry Smith   PetscFunctionReturn(0);
3244b9ad928SBarry Smith }
3254b9ad928SBarry Smith 
3269dcbbd2bSBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
3270d353602SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
3289dcbbd2bSBarry Smith 
3294b9ad928SBarry Smith #undef __FUNCT__
3304b9ad928SBarry Smith #define __FUNCT__ "PCView_MG"
3316849ba73SBarry Smith static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
3324b9ad928SBarry Smith {
3339dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
334dfbe8321SBarry Smith   PetscErrorCode ierr;
33579416396SBarry Smith   PetscInt       levels = mg[0]->levels,i;
33632077d6dSBarry Smith   PetscTruth     iascii;
3374b9ad928SBarry Smith 
3384b9ad928SBarry Smith   PetscFunctionBegin;
33932077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
34032077d6dSBarry Smith   if (iascii) {
341a8e6722dSMatthew Knepley     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg[0]->am],levels,(mg[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr);
342d792f543SMatthew Knepley     if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
343d792f543SMatthew Knepley       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg[0]->cyclesperpcapply);CHKERRQ(ierr);
344d792f543SMatthew Knepley     }
345c2be2410SBarry Smith     if (mg[0]->galerkin) {
346c2be2410SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
347c2be2410SBarry Smith     }
3484b9ad928SBarry Smith     for (i=0; i<levels; i++) {
349b03c7568SBarry Smith       if (!i) {
350a8e6722dSMatthew Knepley         ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D presmooths=%D postsmooths=%D -----\n",i,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr);
351b03c7568SBarry Smith       } else {
352a8e6722dSMatthew Knepley         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg[i]->default_smoothd);CHKERRQ(ierr);
353b03c7568SBarry Smith       }
3544b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3554b9ad928SBarry Smith       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
3564b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
357b03c7568SBarry Smith       if (i && mg[i]->smoothd == mg[i]->smoothu) {
3584b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
359b03c7568SBarry Smith       } else if (i){
360a8e6722dSMatthew Knepley         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg[i]->default_smoothu);CHKERRQ(ierr);
3614b9ad928SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3624b9ad928SBarry Smith         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
3634b9ad928SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3644b9ad928SBarry Smith       }
3654b9ad928SBarry Smith     }
3664b9ad928SBarry Smith   } else {
36779a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
3684b9ad928SBarry Smith   }
3694b9ad928SBarry Smith   PetscFunctionReturn(0);
3704b9ad928SBarry Smith }
3714b9ad928SBarry Smith 
3724b9ad928SBarry Smith /*
3734b9ad928SBarry Smith     Calls setup for the KSP on each level
3744b9ad928SBarry Smith */
3754b9ad928SBarry Smith #undef __FUNCT__
3764b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_MG"
3776849ba73SBarry Smith static PetscErrorCode PCSetUp_MG(PC pc)
3784b9ad928SBarry Smith {
3799dcbbd2bSBarry Smith   PC_MG                   **mg = (PC_MG**)pc->data;
380dfbe8321SBarry Smith   PetscErrorCode          ierr;
38179416396SBarry Smith   PetscInt                i,n = mg[0]->levels;
38277122347SBarry Smith   PC                      cpc,mpc;
383906ed7ccSBarry Smith   PetscTruth              preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset;
38423d894e5SBarry Smith   PetscViewerASCIIMonitor ascii;
38523d894e5SBarry Smith   PetscViewer             viewer = PETSC_NULL;
3864b9ad928SBarry Smith   MPI_Comm                comm;
387c2be2410SBarry Smith   Mat                     dA,dB;
388c2be2410SBarry Smith   MatStructure            uflag;
3890a6bb862SBarry Smith   Vec                     tvec;
3904b9ad928SBarry Smith 
3914b9ad928SBarry Smith   PetscFunctionBegin;
392852665d3SBarry Smith 
39343fb2f97SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
39443fb2f97SBarry Smith   /* so use those from global PC */
39543fb2f97SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
396906ed7ccSBarry Smith   ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
39743fb2f97SBarry Smith   ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
39877122347SBarry Smith   ierr = KSPGetPC(mg[n-1]->smoothd,&mpc);CHKERRQ(ierr);
39977122347SBarry Smith   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2))) {
400852665d3SBarry 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);
4011cfe3bddSBarry Smith     ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
402852665d3SBarry Smith   }
403852665d3SBarry Smith 
404852665d3SBarry Smith   if (mg[0]->galerkin) {
405852665d3SBarry Smith     Mat B;
406852665d3SBarry Smith     mg[0]->galerkinused = PETSC_TRUE;
407852665d3SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
408852665d3SBarry Smith     ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
409852665d3SBarry Smith     if (!pc->setupcalled) {
410852665d3SBarry Smith       for (i=n-2; i>-1; i--) {
411852665d3SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
412852665d3SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
413906ed7ccSBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
414852665d3SBarry Smith         dB   = B;
415852665d3SBarry Smith       }
416906ed7ccSBarry Smith       ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);
417852665d3SBarry Smith     } else {
418852665d3SBarry Smith       for (i=n-2; i>-1; i--) {
419906ed7ccSBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
420852665d3SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
421852665d3SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
422852665d3SBarry Smith         dB   = B;
423852665d3SBarry Smith       }
424852665d3SBarry Smith     }
425852665d3SBarry Smith   }
426852665d3SBarry Smith 
427958c9bccSBarry Smith   if (!pc->setupcalled) {
4284b9ad928SBarry Smith     ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr);
4294b9ad928SBarry Smith 
430b03c7568SBarry Smith     for (i=0; i<n; i++) {
4314b9ad928SBarry Smith       if (monitor) {
4324b9ad928SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
43323d894e5SBarry Smith         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
43423d894e5SBarry Smith         ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
4354b9ad928SBarry Smith       }
4364b9ad928SBarry Smith       ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
4374b9ad928SBarry Smith     }
438b03c7568SBarry Smith     for (i=1; i<n; i++) {
439a98fc643SBarry Smith       if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) {
4404b9ad928SBarry Smith         if (monitor) {
4414b9ad928SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
44223d894e5SBarry Smith           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
44323d894e5SBarry Smith           ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
4444b9ad928SBarry Smith         }
4454b9ad928SBarry Smith         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
4464b9ad928SBarry Smith       }
4474b9ad928SBarry Smith     }
448fccaa45eSBarry Smith     for (i=1; i<n; i++) {
4490cace4b0SBarry Smith       if (!mg[i]->residual) {
4500cace4b0SBarry Smith         Mat mat;
4510cace4b0SBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
4520cace4b0SBarry Smith         ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
4530cace4b0SBarry Smith       }
454fccaa45eSBarry Smith       if (mg[i]->restrct && !mg[i]->interpolate) {
455aea2a34eSBarry Smith         ierr = PCMGSetInterpolation(pc,i,mg[i]->restrct);CHKERRQ(ierr);
456fccaa45eSBarry Smith       }
457fccaa45eSBarry Smith       if (!mg[i]->restrct && mg[i]->interpolate) {
45897177400SBarry Smith         ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr);
459fccaa45eSBarry Smith       }
460fccaa45eSBarry Smith #if defined(PETSC_USE_DEBUG)
461fccaa45eSBarry Smith       if (!mg[i]->restrct || !mg[i]->interpolate) {
462fccaa45eSBarry Smith         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
463fccaa45eSBarry Smith       }
464fccaa45eSBarry Smith #endif
465fccaa45eSBarry Smith     }
4660a6bb862SBarry Smith     for (i=0; i<n-1; i++) {
46737b0e6c0SBarry Smith       if (!mg[i]->b) {
468906ed7ccSBarry Smith         Vec *vec;
469906ed7ccSBarry Smith         ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
470906ed7ccSBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
4714fdc791cSSatish Balay         ierr = VecDestroy(*vec);CHKERRQ(ierr);
472906ed7ccSBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
47337b0e6c0SBarry Smith       }
4746ca4d86aSHong Zhang       if (!mg[i]->r && i) {
4750a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
47697177400SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
4770a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4780a6bb862SBarry Smith       }
4790a6bb862SBarry Smith       if (!mg[i]->x) {
4800a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
48197177400SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
4820a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4830a6bb862SBarry Smith       }
4840a6bb862SBarry Smith     }
485dfef70bdSHong Zhang     if (n != 1 && !mg[n-1]->r) {
486dfef70bdSHong Zhang       /* PCMGSetR() on the finest level if user did not supply it */
4870f77fa87SBarry Smith       Vec *vec;
4880f77fa87SBarry Smith       ierr = KSPGetVecs(mg[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
4890f77fa87SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
490797d2ccbSBarry Smith       ierr = VecDestroy(*vec);CHKERRQ(ierr);
4910f77fa87SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
4920f77fa87SBarry Smith     }
4934b9ad928SBarry Smith   }
4944b9ad928SBarry Smith 
495c2be2410SBarry Smith 
4964b9ad928SBarry Smith   for (i=1; i<n; i++) {
497b03c7568SBarry Smith     if (mg[i]->smoothu == mg[i]->smoothd) {
498b03c7568SBarry Smith       /* if doing only down then initial guess is zero */
4994b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
500b03c7568SBarry Smith     }
50132cf1786SBarry Smith     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5024b9ad928SBarry Smith     ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
50332cf1786SBarry Smith     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5044b9ad928SBarry Smith   }
505b03c7568SBarry Smith   for (i=1; i<n; i++) {
5064b9ad928SBarry Smith     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
507906ed7ccSBarry Smith       Mat          downmat,downpmat;
50897f1f81fSBarry Smith       MatStructure matflag;
509906ed7ccSBarry Smith       PetscTruth   opsset;
51097f1f81fSBarry Smith 
51197f1f81fSBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
512906ed7ccSBarry Smith       ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
513906ed7ccSBarry Smith       if (!opsset) {
514906ed7ccSBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
51597f1f81fSBarry Smith         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
51697f1f81fSBarry Smith       }
51797f1f81fSBarry Smith 
5184b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
51932cf1786SBarry Smith       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5204b9ad928SBarry Smith       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
52132cf1786SBarry Smith       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5224b9ad928SBarry Smith     }
5234b9ad928SBarry Smith   }
5244b9ad928SBarry Smith 
5254b9ad928SBarry Smith   /*
5264b9ad928SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
5274b9ad928SBarry Smith   */
5284b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
5294b9ad928SBarry Smith   if (preonly) {
5304b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
5314b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
53268eff7e6SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
53368eff7e6SBarry Smith     if (!lu && !redundant && !cholesky) {
5344b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
5354b9ad928SBarry Smith     }
5364b9ad928SBarry Smith   }
5374b9ad928SBarry Smith 
538958c9bccSBarry Smith   if (!pc->setupcalled) {
5394b9ad928SBarry Smith     if (monitor) {
5404b9ad928SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
54123d894e5SBarry Smith       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
54223d894e5SBarry Smith       ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
5434b9ad928SBarry Smith     }
5444b9ad928SBarry Smith     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
5454b9ad928SBarry Smith   }
5464b9ad928SBarry Smith 
54732cf1786SBarry Smith   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5484b9ad928SBarry Smith   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
54932cf1786SBarry Smith   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5504b9ad928SBarry Smith 
5514b9ad928SBarry Smith   /*
5526805f65bSBarry Smith      Dump the interpolation/restriction matrices plus the
5534b9ad928SBarry Smith    Jacobian/stiffness on each level. This allows Matlab users to
5546805f65bSBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
5556805f65bSBarry Smith 
5566805f65bSBarry Smith    Only support one or the other at the same time.
5576805f65bSBarry Smith   */
5586805f65bSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
5597adad957SLisandro Dalcin   ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr);
5604b9ad928SBarry Smith   if (dump) {
5617adad957SLisandro Dalcin     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
5624b9ad928SBarry Smith   }
563c45a1595SBarry Smith #endif
5647adad957SLisandro Dalcin   ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr);
565c2be2410SBarry Smith   if (dump) {
5667adad957SLisandro Dalcin     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
5676805f65bSBarry Smith   }
5686805f65bSBarry Smith 
5696805f65bSBarry Smith   if (viewer) {
570c2be2410SBarry Smith     for (i=1; i<n; i++) {
5716805f65bSBarry Smith       ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr);
572c2be2410SBarry Smith     }
573c2be2410SBarry Smith     for (i=0; i<n; i++) {
574c2be2410SBarry Smith       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
5756805f65bSBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
576c2be2410SBarry Smith     }
577c2be2410SBarry Smith   }
5784b9ad928SBarry Smith   PetscFunctionReturn(0);
5794b9ad928SBarry Smith }
5804b9ad928SBarry Smith 
5814b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
5824b9ad928SBarry Smith 
5834b9ad928SBarry Smith #undef __FUNCT__
5849dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
5854b9ad928SBarry Smith /*@C
58697177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
5874b9ad928SBarry Smith    Must be called before any other MG routine.
5884b9ad928SBarry Smith 
5894b9ad928SBarry Smith    Collective on PC
5904b9ad928SBarry Smith 
5914b9ad928SBarry Smith    Input Parameters:
5924b9ad928SBarry Smith +  pc - the preconditioner context
5934b9ad928SBarry Smith .  levels - the number of levels
5944b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
5954b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
5964b9ad928SBarry Smith 
5974b9ad928SBarry Smith    Level: intermediate
5984b9ad928SBarry Smith 
5994b9ad928SBarry Smith    Notes:
6004b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
6014b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
6024b9ad928SBarry Smith 
6034b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
6044b9ad928SBarry Smith 
60597177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
6064b9ad928SBarry Smith @*/
60797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
6084b9ad928SBarry Smith {
609dfbe8321SBarry Smith   PetscErrorCode ierr;
610ada7143aSSatish Balay   PC_MG          **mg=0;
6114b9ad928SBarry Smith 
6124b9ad928SBarry Smith   PetscFunctionBegin;
6134482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6144b9ad928SBarry Smith 
6154b9ad928SBarry Smith   if (pc->data) {
6161302d50aSBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
61797177400SBarry Smith     make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
6184b9ad928SBarry Smith   }
6197adad957SLisandro Dalcin   ierr                     = PCMGCreate_Private(((PetscObject)pc)->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
6209dcbbd2bSBarry Smith   mg[0]->am                = PC_MG_MULTIPLICATIVE;
6214b9ad928SBarry Smith   pc->data                 = (void*)mg;
6224b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
6234b9ad928SBarry Smith   PetscFunctionReturn(0);
6244b9ad928SBarry Smith }
6254b9ad928SBarry Smith 
6264b9ad928SBarry Smith #undef __FUNCT__
6279dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
6284b9ad928SBarry Smith /*@
62997177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
6304b9ad928SBarry Smith 
6314b9ad928SBarry Smith    Not Collective
6324b9ad928SBarry Smith 
6334b9ad928SBarry Smith    Input Parameter:
6344b9ad928SBarry Smith .  pc - the preconditioner context
6354b9ad928SBarry Smith 
6364b9ad928SBarry Smith    Output parameter:
6374b9ad928SBarry Smith .  levels - the number of levels
6384b9ad928SBarry Smith 
6394b9ad928SBarry Smith    Level: advanced
6404b9ad928SBarry Smith 
6414b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
6424b9ad928SBarry Smith 
64397177400SBarry Smith .seealso: PCMGSetLevels()
6444b9ad928SBarry Smith @*/
64597177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
6464b9ad928SBarry Smith {
6479dcbbd2bSBarry Smith   PC_MG  **mg;
6484b9ad928SBarry Smith 
6494b9ad928SBarry Smith   PetscFunctionBegin;
6504482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6514482741eSBarry Smith   PetscValidIntPointer(levels,2);
6524b9ad928SBarry Smith 
6539dcbbd2bSBarry Smith   mg      = (PC_MG**)pc->data;
6544b9ad928SBarry Smith   *levels = mg[0]->levels;
6554b9ad928SBarry Smith   PetscFunctionReturn(0);
6564b9ad928SBarry Smith }
6574b9ad928SBarry Smith 
6584b9ad928SBarry Smith #undef __FUNCT__
6599dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
6604b9ad928SBarry Smith /*@
66197177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
6624b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
6634b9ad928SBarry Smith 
6644b9ad928SBarry Smith    Collective on PC
6654b9ad928SBarry Smith 
6664b9ad928SBarry Smith    Input Parameters:
6674b9ad928SBarry Smith +  pc - the preconditioner context
6689dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
6699dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
6704b9ad928SBarry Smith 
6714b9ad928SBarry Smith    Options Database Key:
6724b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
6734b9ad928SBarry Smith    additive, full, kaskade
6744b9ad928SBarry Smith 
6754b9ad928SBarry Smith    Level: advanced
6764b9ad928SBarry Smith 
6774b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
6784b9ad928SBarry Smith 
67997177400SBarry Smith .seealso: PCMGSetLevels()
6804b9ad928SBarry Smith @*/
6819dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
6824b9ad928SBarry Smith {
6839dcbbd2bSBarry Smith   PC_MG **mg;
6844b9ad928SBarry Smith 
6854b9ad928SBarry Smith   PetscFunctionBegin;
6864482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6879dcbbd2bSBarry Smith   mg = (PC_MG**)pc->data;
6884b9ad928SBarry Smith 
6894b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
6904b9ad928SBarry Smith   mg[0]->am = form;
6919dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
6924b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
6934b9ad928SBarry Smith   PetscFunctionReturn(0);
6944b9ad928SBarry Smith }
6954b9ad928SBarry Smith 
6964b9ad928SBarry Smith #undef __FUNCT__
6970d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
6984b9ad928SBarry Smith /*@
6990d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
7004b9ad928SBarry Smith    complicated cycling.
7014b9ad928SBarry Smith 
7024b9ad928SBarry Smith    Collective on PC
7034b9ad928SBarry Smith 
7044b9ad928SBarry Smith    Input Parameters:
705c2be2410SBarry Smith +  pc - the multigrid context
7060d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
7074b9ad928SBarry Smith 
7084b9ad928SBarry Smith    Options Database Key:
7090d353602SBarry Smith $  -pc_mg_cycle_type v or w
7104b9ad928SBarry Smith 
7114b9ad928SBarry Smith    Level: advanced
7124b9ad928SBarry Smith 
7134b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7144b9ad928SBarry Smith 
7150d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
7164b9ad928SBarry Smith @*/
7170d353602SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n)
7184b9ad928SBarry Smith {
7199dcbbd2bSBarry Smith   PC_MG    **mg;
72079416396SBarry Smith   PetscInt i,levels;
7214b9ad928SBarry Smith 
7224b9ad928SBarry Smith   PetscFunctionBegin;
7234482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7249dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
7254b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7264b9ad928SBarry Smith   levels = mg[0]->levels;
7274b9ad928SBarry Smith 
7284b9ad928SBarry Smith   for (i=0; i<levels; i++) {
7294b9ad928SBarry Smith     mg[i]->cycles  = n;
7304b9ad928SBarry Smith   }
7314b9ad928SBarry Smith   PetscFunctionReturn(0);
7324b9ad928SBarry Smith }
7334b9ad928SBarry Smith 
7344b9ad928SBarry Smith #undef __FUNCT__
7358cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
7368cc2d5dfSBarry Smith /*@
7378cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
7388cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
7398cc2d5dfSBarry Smith 
7408cc2d5dfSBarry Smith    Collective on PC
7418cc2d5dfSBarry Smith 
7428cc2d5dfSBarry Smith    Input Parameters:
7438cc2d5dfSBarry Smith +  pc - the multigrid context
7448cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
7458cc2d5dfSBarry Smith 
7468cc2d5dfSBarry Smith    Options Database Key:
7478cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
7488cc2d5dfSBarry Smith 
7498cc2d5dfSBarry Smith    Level: advanced
7508cc2d5dfSBarry Smith 
7518cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
7528cc2d5dfSBarry Smith 
7538cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7548cc2d5dfSBarry Smith 
7558cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
7568cc2d5dfSBarry Smith @*/
7578cc2d5dfSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
7588cc2d5dfSBarry Smith {
7598cc2d5dfSBarry Smith   PC_MG    **mg;
7608cc2d5dfSBarry Smith   PetscInt i,levels;
7618cc2d5dfSBarry Smith 
7628cc2d5dfSBarry Smith   PetscFunctionBegin;
7638cc2d5dfSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7648cc2d5dfSBarry Smith   mg     = (PC_MG**)pc->data;
7658cc2d5dfSBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7668cc2d5dfSBarry Smith   levels = mg[0]->levels;
7678cc2d5dfSBarry Smith 
7688cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
7698cc2d5dfSBarry Smith     mg[i]->cyclesperpcapply  = n;
7708cc2d5dfSBarry Smith   }
7718cc2d5dfSBarry Smith   PetscFunctionReturn(0);
7728cc2d5dfSBarry Smith }
7738cc2d5dfSBarry Smith 
7748cc2d5dfSBarry Smith #undef __FUNCT__
7759dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
776c2be2410SBarry Smith /*@
77797177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
778c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
779c2be2410SBarry Smith 
780c2be2410SBarry Smith    Collective on PC
781c2be2410SBarry Smith 
782c2be2410SBarry Smith    Input Parameters:
7833fc8bf9cSBarry Smith .  pc - the multigrid context
784c2be2410SBarry Smith 
785c2be2410SBarry Smith    Options Database Key:
786c2be2410SBarry Smith $  -pc_mg_galerkin
787c2be2410SBarry Smith 
788c2be2410SBarry Smith    Level: intermediate
789c2be2410SBarry Smith 
790c2be2410SBarry Smith .keywords: MG, set, Galerkin
791c2be2410SBarry Smith 
7923fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
7933fc8bf9cSBarry Smith 
794c2be2410SBarry Smith @*/
79597177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
796c2be2410SBarry Smith {
7979dcbbd2bSBarry Smith   PC_MG    **mg;
798c2be2410SBarry Smith   PetscInt i,levels;
799c2be2410SBarry Smith 
800c2be2410SBarry Smith   PetscFunctionBegin;
801c2be2410SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8029dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
803c2be2410SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
804c2be2410SBarry Smith   levels = mg[0]->levels;
805c2be2410SBarry Smith 
806c2be2410SBarry Smith   for (i=0; i<levels; i++) {
807c2be2410SBarry Smith     mg[i]->galerkin = PETSC_TRUE;
808c2be2410SBarry Smith   }
809c2be2410SBarry Smith   PetscFunctionReturn(0);
810c2be2410SBarry Smith }
811c2be2410SBarry Smith 
812c2be2410SBarry Smith #undef __FUNCT__
8133fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
8143fc8bf9cSBarry Smith /*@
8153fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
8163fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
8173fc8bf9cSBarry Smith 
8183fc8bf9cSBarry Smith    Not Collective
8193fc8bf9cSBarry Smith 
8203fc8bf9cSBarry Smith    Input Parameter:
8213fc8bf9cSBarry Smith .  pc - the multigrid context
8223fc8bf9cSBarry Smith 
8233fc8bf9cSBarry Smith    Output Parameter:
8243fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
8253fc8bf9cSBarry Smith 
8263fc8bf9cSBarry Smith    Options Database Key:
8273fc8bf9cSBarry Smith $  -pc_mg_galerkin
8283fc8bf9cSBarry Smith 
8293fc8bf9cSBarry Smith    Level: intermediate
8303fc8bf9cSBarry Smith 
8313fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
8323fc8bf9cSBarry Smith 
8333fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
8343fc8bf9cSBarry Smith 
8353fc8bf9cSBarry Smith @*/
8363fc8bf9cSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
8373fc8bf9cSBarry Smith {
8383fc8bf9cSBarry Smith   PC_MG    **mg;
8393fc8bf9cSBarry Smith 
8403fc8bf9cSBarry Smith   PetscFunctionBegin;
8413fc8bf9cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8423fc8bf9cSBarry Smith   mg     = (PC_MG**)pc->data;
8433fc8bf9cSBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
8443fc8bf9cSBarry Smith   *galerkin = mg[0]->galerkin;
8453fc8bf9cSBarry Smith   PetscFunctionReturn(0);
8463fc8bf9cSBarry Smith }
8473fc8bf9cSBarry Smith 
8483fc8bf9cSBarry Smith #undef __FUNCT__
8499dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
8504b9ad928SBarry Smith /*@
85197177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
85297177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
8534b9ad928SBarry Smith    pre-smoothing steps on different levels.
8544b9ad928SBarry Smith 
8554b9ad928SBarry Smith    Collective on PC
8564b9ad928SBarry Smith 
8574b9ad928SBarry Smith    Input Parameters:
8584b9ad928SBarry Smith +  mg - the multigrid context
8594b9ad928SBarry Smith -  n - the number of smoothing steps
8604b9ad928SBarry Smith 
8614b9ad928SBarry Smith    Options Database Key:
8624b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
8634b9ad928SBarry Smith 
8644b9ad928SBarry Smith    Level: advanced
8654b9ad928SBarry Smith 
8664b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
8674b9ad928SBarry Smith 
86897177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
8694b9ad928SBarry Smith @*/
87097177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
8714b9ad928SBarry Smith {
8729dcbbd2bSBarry Smith   PC_MG          **mg;
8736849ba73SBarry Smith   PetscErrorCode ierr;
87479416396SBarry Smith   PetscInt       i,levels;
8754b9ad928SBarry Smith 
8764b9ad928SBarry Smith   PetscFunctionBegin;
8774482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8789dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
8794b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
8804b9ad928SBarry Smith   levels = mg[0]->levels;
8814b9ad928SBarry Smith 
882b05257ddSBarry Smith   for (i=1; i<levels; i++) {
8834b9ad928SBarry Smith     /* make sure smoother up and down are different */
88497177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
8854b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
8864b9ad928SBarry Smith     mg[i]->default_smoothd = n;
8874b9ad928SBarry Smith   }
8884b9ad928SBarry Smith   PetscFunctionReturn(0);
8894b9ad928SBarry Smith }
8904b9ad928SBarry Smith 
8914b9ad928SBarry Smith #undef __FUNCT__
8929dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
8934b9ad928SBarry Smith /*@
89497177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
89597177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
8964b9ad928SBarry Smith    post-smoothing steps on different levels.
8974b9ad928SBarry Smith 
8984b9ad928SBarry Smith    Collective on PC
8994b9ad928SBarry Smith 
9004b9ad928SBarry Smith    Input Parameters:
9014b9ad928SBarry Smith +  mg - the multigrid context
9024b9ad928SBarry Smith -  n - the number of smoothing steps
9034b9ad928SBarry Smith 
9044b9ad928SBarry Smith    Options Database Key:
9054b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
9064b9ad928SBarry Smith 
9074b9ad928SBarry Smith    Level: advanced
9084b9ad928SBarry Smith 
9094b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
910a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
9114b9ad928SBarry Smith 
9124b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
9134b9ad928SBarry Smith 
91497177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
9154b9ad928SBarry Smith @*/
91697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
9174b9ad928SBarry Smith {
9189dcbbd2bSBarry Smith   PC_MG          **mg;
9196849ba73SBarry Smith   PetscErrorCode ierr;
92079416396SBarry Smith   PetscInt       i,levels;
9214b9ad928SBarry Smith 
9224b9ad928SBarry Smith   PetscFunctionBegin;
9234482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
9249dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
9254b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
9264b9ad928SBarry Smith   levels = mg[0]->levels;
9274b9ad928SBarry Smith 
9284b9ad928SBarry Smith   for (i=1; i<levels; i++) {
9294b9ad928SBarry Smith     /* make sure smoother up and down are different */
93097177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
9314b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
9324b9ad928SBarry Smith     mg[i]->default_smoothu = n;
9334b9ad928SBarry Smith   }
9344b9ad928SBarry Smith   PetscFunctionReturn(0);
9354b9ad928SBarry Smith }
9364b9ad928SBarry Smith 
9374b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
9384b9ad928SBarry Smith 
9393b09bd56SBarry Smith /*MC
940ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
9413b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
9423b09bd56SBarry Smith 
9433b09bd56SBarry Smith    Options Database Keys:
9443b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
9450d353602SBarry Smith .  -pc_mg_cycles v or w
94679416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
9473b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
9483b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
9493b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
9503b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
95168eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
9523b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
9533b09bd56SBarry Smith                         to the Socket viewer for reading from Matlab.
9543b09bd56SBarry Smith 
9553b09bd56SBarry Smith    Notes:
9563b09bd56SBarry Smith 
9573b09bd56SBarry Smith    Level: intermediate
9583b09bd56SBarry Smith 
9598f87f92bSBarry Smith    Concepts: multigrid/multilevel
9603b09bd56SBarry Smith 
9613b09bd56SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
9620d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
96397177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
96497177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
9650d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
9663b09bd56SBarry Smith M*/
9673b09bd56SBarry Smith 
9684b9ad928SBarry Smith EXTERN_C_BEGIN
9694b9ad928SBarry Smith #undef __FUNCT__
9704b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
971dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
9724b9ad928SBarry Smith {
9734b9ad928SBarry Smith   PetscFunctionBegin;
9744b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
9754b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
9764b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
9774b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
9784b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
9794b9ad928SBarry Smith 
9804b9ad928SBarry Smith   pc->data                = (void*)0;
9814b9ad928SBarry Smith   PetscFunctionReturn(0);
9824b9ad928SBarry Smith }
9834b9ad928SBarry Smith EXTERN_C_END
984