xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 90d69ab76acc68be38016f1844fbdbca1424683c)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
54b9ad928SBarry Smith */
67c4f633dSBarry 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"
114d0a8057SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG **mglevels,PCRichardsonConvergedReason *reason)
124b9ad928SBarry Smith {
139dcbbd2bSBarry Smith   PC_MG          *mg = *mglevels,*mgc;
146849ba73SBarry Smith   PetscErrorCode ierr;
15e0ef7d32SBarry Smith   PetscInt       cycles = (mg->level == 1) ? 1 : (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 */
284d0a8057SBarry 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) {
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 
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--) {
494d0a8057SBarry 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);
961cee3971SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mg[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
9779416396SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr);
984b9ad928SBarry Smith     ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr);
994b9ad928SBarry Smith 
1004b9ad928SBarry Smith     /* do special stuff for coarse grid */
1014b9ad928SBarry Smith     if (!i && levels > 1) {
1024b9ad928SBarry Smith       ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
1034b9ad928SBarry Smith 
1044b9ad928SBarry Smith       /* coarse solve is (redundant) LU by default */
1054b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
1064b9ad928SBarry Smith       ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr);
1074b9ad928SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1084b9ad928SBarry Smith       if (size > 1) {
1094b9ad928SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
1100d7810c8SBarry Smith       } else {
1114b9ad928SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
1120d7810c8SBarry Smith       }
1134b9ad928SBarry Smith 
1144b9ad928SBarry Smith     } else {
1152e70eadcSBarry Smith       char tprefix[128];
11613f74950SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
1172e70eadcSBarry Smith       ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr);
1184b9ad928SBarry Smith     }
11952e6d16bSBarry Smith     ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr);
1204b9ad928SBarry Smith     mg[i]->smoothu             = mg[i]->smoothd;
1214b9ad928SBarry Smith     mg[i]->rtol                = 0.0;
12270441072SBarry Smith     mg[i]->abstol              = 0.0;
1234b9ad928SBarry Smith     mg[i]->dtol                = 0.0;
1244b9ad928SBarry Smith     mg[i]->ttol                = 0.0;
12532cf1786SBarry Smith     mg[i]->eventsmoothsetup    = 0;
12632cf1786SBarry Smith     mg[i]->eventsmoothsolve    = 0;
12732cf1786SBarry Smith     mg[i]->eventresidual       = 0;
12832cf1786SBarry Smith     mg[i]->eventinterprestrict = 0;
1298cc2d5dfSBarry Smith     mg[i]->cyclesperpcapply    = 1;
1304b9ad928SBarry Smith   }
1314b9ad928SBarry Smith   *result = mg;
1324b9ad928SBarry Smith   PetscFunctionReturn(0);
1334b9ad928SBarry Smith }
1344b9ad928SBarry Smith 
1354b9ad928SBarry Smith #undef __FUNCT__
1364b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_MG"
1377233f9f0SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
1384b9ad928SBarry Smith {
1399dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1406849ba73SBarry Smith   PetscErrorCode ierr;
14138f2d2fdSLisandro Dalcin   PetscInt       i,n;
1424b9ad928SBarry Smith 
1434b9ad928SBarry Smith   PetscFunctionBegin;
14438f2d2fdSLisandro Dalcin   if (!mg) PetscFunctionReturn(0);
14514ca34e6SBarry Smith 
14638f2d2fdSLisandro Dalcin   n = mg[0]->levels;
147fccaa45eSBarry Smith   for (i=0; i<n-1; i++) {
148630ba2eeSBarry Smith     if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);}
149fccaa45eSBarry Smith     if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);}
150fccaa45eSBarry Smith     if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);}
151fccaa45eSBarry Smith     if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);}
152fccaa45eSBarry Smith     if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);}
153fccaa45eSBarry Smith   }
154fccaa45eSBarry Smith 
1554b9ad928SBarry Smith   for (i=0; i<n; i++) {
1564b9ad928SBarry Smith     if (mg[i]->smoothd != mg[i]->smoothu) {
1574b9ad928SBarry Smith       ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr);
1584b9ad928SBarry Smith     }
1594b9ad928SBarry Smith     ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr);
1604b9ad928SBarry Smith     ierr = PetscFree(mg[i]);CHKERRQ(ierr);
1614b9ad928SBarry Smith   }
1624b9ad928SBarry Smith   ierr = PetscFree(mg);CHKERRQ(ierr);
1634b9ad928SBarry Smith   PetscFunctionReturn(0);
1644b9ad928SBarry Smith }
1654b9ad928SBarry Smith 
1664b9ad928SBarry Smith 
1674b9ad928SBarry Smith 
1689dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**);
1691e2582c4SBarry Smith EXTERN PetscErrorCode PCMGFCycle_Private(PC,PC_MG**);
1709dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**);
1714b9ad928SBarry Smith 
1724b9ad928SBarry Smith /*
1734b9ad928SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
1744b9ad928SBarry Smith              or full cycle of multigrid.
1754b9ad928SBarry Smith 
1764b9ad928SBarry Smith   Note:
1779dcbbd2bSBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
1784b9ad928SBarry Smith */
1794b9ad928SBarry Smith #undef __FUNCT__
1804b9ad928SBarry Smith #define __FUNCT__ "PCApply_MG"
1816849ba73SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
1824b9ad928SBarry Smith {
1839dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
1846849ba73SBarry Smith   PetscErrorCode ierr;
1858cc2d5dfSBarry Smith   PetscInt       levels = mg[0]->levels,i;
1864b9ad928SBarry Smith 
1874b9ad928SBarry Smith   PetscFunctionBegin;
1884b9ad928SBarry Smith   mg[levels-1]->b = b;
1894b9ad928SBarry Smith   mg[levels-1]->x = x;
1909dcbbd2bSBarry Smith   if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
191efb30889SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
1928cc2d5dfSBarry Smith     for (i=0; i<mg[0]->cyclesperpcapply; i++) {
1931e2582c4SBarry Smith       ierr = PCMGMCycle_Private(pc,mg+levels-1,PETSC_NULL);CHKERRQ(ierr);
1944b9ad928SBarry Smith     }
1958cc2d5dfSBarry Smith   }
1969dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_ADDITIVE) {
1979dcbbd2bSBarry Smith     ierr = PCMGACycle_Private(mg);CHKERRQ(ierr);
1984b9ad928SBarry Smith   }
1999dcbbd2bSBarry Smith   else if (mg[0]->am == PC_MG_KASKADE) {
2009dcbbd2bSBarry Smith     ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr);
2014b9ad928SBarry Smith   }
2024b9ad928SBarry Smith   else {
2031e2582c4SBarry Smith     ierr = PCMGFCycle_Private(pc,mg);CHKERRQ(ierr);
2044b9ad928SBarry Smith   }
2054b9ad928SBarry Smith   PetscFunctionReturn(0);
2064b9ad928SBarry Smith }
2074b9ad928SBarry Smith 
2084b9ad928SBarry Smith #undef __FUNCT__
2094b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
2104d0a8057SBarry 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)
2114b9ad928SBarry Smith {
2129dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
213dfbe8321SBarry Smith   PetscErrorCode ierr;
2144d0a8057SBarry Smith   PetscInt       levels = mg[0]->levels,i;
2154b9ad928SBarry Smith 
2164b9ad928SBarry Smith   PetscFunctionBegin;
2174b9ad928SBarry Smith   mg[levels-1]->b    = b;
2184b9ad928SBarry Smith   mg[levels-1]->x    = x;
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith   mg[levels-1]->rtol = rtol;
22170441072SBarry Smith   mg[levels-1]->abstol = abstol;
2224b9ad928SBarry Smith   mg[levels-1]->dtol = dtol;
2234b9ad928SBarry Smith   if (rtol) {
2244b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
2254b9ad928SBarry Smith     PetscReal rnorm;
2264b9ad928SBarry Smith     ierr               = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr);
2274b9ad928SBarry Smith     ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
22870441072SBarry Smith     mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol);
22970441072SBarry Smith   } else if (abstol) {
23070441072SBarry Smith     mg[levels-1]->ttol = abstol;
2314b9ad928SBarry Smith   } else {
2324b9ad928SBarry Smith     mg[levels-1]->ttol = 0.0;
2334b9ad928SBarry Smith   }
2344b9ad928SBarry Smith 
2354d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
2364d0a8057SBarry Smith      stop prematurely do to small residual */
2374d0a8057SBarry Smith   for (i=1; i<levels; i++) {
2384d0a8057SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2394d0a8057SBarry Smith     if (mg[i]->smoothu != mg[i]->smoothd) {
2404d0a8057SBarry Smith       ierr = KSPSetTolerances(mg[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
2414b9ad928SBarry Smith     }
2424d0a8057SBarry Smith   }
2434d0a8057SBarry Smith 
2444d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
2454d0a8057SBarry Smith   for (i=0; i<its; i++) {
2464d0a8057SBarry Smith     ierr = PCMGMCycle_Private(pc,mg+levels-1,reason);CHKERRQ(ierr);
2474d0a8057SBarry Smith     if (*reason) break;
2484d0a8057SBarry Smith   }
2494d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
2504d0a8057SBarry Smith   *outits = i;
2514b9ad928SBarry Smith   PetscFunctionReturn(0);
2524b9ad928SBarry Smith }
2534b9ad928SBarry Smith 
2544b9ad928SBarry Smith #undef __FUNCT__
2554b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
2566ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_MG(PC pc)
2574b9ad928SBarry Smith {
258dfbe8321SBarry Smith   PetscErrorCode ierr;
2598cc2d5dfSBarry Smith   PetscInt       m,levels = 1,cycles;
2604b9ad928SBarry Smith   PetscTruth     flg;
2619dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
26290da95b0SMatthew Knepley   PCMGType       mgtype = PC_MG_ADDITIVE;
2631aa34eeaSBarry Smith   PCMGCycleType  mgctype;
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith   PetscFunctionBegin;
2664b9ad928SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
2674b9ad928SBarry Smith     if (!pc->data) {
2689dcbbd2bSBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
26997177400SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
270cf502942SBarry Smith       mg = (PC_MG**)pc->data;
2714b9ad928SBarry Smith     }
272f2070a76SMatthew Knepley     mgctype = (PCMGCycleType) mg[0]->cycles;
2730d353602SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
2744b9ad928SBarry Smith     if (flg) {
2750d353602SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
2760d353602SBarry Smith     };
277*90d69ab7SBarry Smith     flg  = PETSC_FALSE;
278*90d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
279c2be2410SBarry Smith     if (flg) {
28097177400SBarry Smith       ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr);
281c2be2410SBarry Smith     }
2829dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
2834b9ad928SBarry Smith     if (flg) {
28497177400SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
2854b9ad928SBarry Smith     }
2869dcbbd2bSBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
2874b9ad928SBarry Smith     if (flg) {
28897177400SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
2894b9ad928SBarry Smith     }
2909dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
2918cc2d5dfSBarry Smith     if (flg) {
2928cc2d5dfSBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2938cc2d5dfSBarry Smith     }
2948cc2d5dfSBarry Smith     if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
2958cc2d5dfSBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
2968cc2d5dfSBarry Smith       if (flg) {
2978cc2d5dfSBarry Smith 	ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
2988cc2d5dfSBarry Smith       }
2998cc2d5dfSBarry Smith     }
300*90d69ab7SBarry Smith     flg  = PETSC_FALSE;
301*90d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
3024b9ad928SBarry Smith     if (flg) {
3034f5ab15aSBarry Smith       PetscInt i;
3044b9ad928SBarry Smith       char     eventname[128];
3054b9ad928SBarry Smith       if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
3064b9ad928SBarry Smith       levels = mg[0]->levels;
3074b9ad928SBarry Smith       for (i=0; i<levels; i++) {
30832cf1786SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
3098cbcd9ccSBarry Smith         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventsmoothsetup);CHKERRQ(ierr);
31032cf1786SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
3118cbcd9ccSBarry Smith         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventsmoothsolve);CHKERRQ(ierr);
31232cf1786SBarry Smith         if (i) {
31332cf1786SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
314a3bc4eb9SBarry Smith           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventresidual);CHKERRQ(ierr);
31532cf1786SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
3168cbcd9ccSBarry Smith           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventinterprestrict);CHKERRQ(ierr);
31732cf1786SBarry Smith         }
3184b9ad928SBarry Smith       }
3194b9ad928SBarry Smith     }
3204b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
3214b9ad928SBarry Smith   PetscFunctionReturn(0);
3224b9ad928SBarry Smith }
3234b9ad928SBarry Smith 
3249dcbbd2bSBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
3250d353602SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
3269dcbbd2bSBarry Smith 
3274b9ad928SBarry Smith #undef __FUNCT__
3284b9ad928SBarry Smith #define __FUNCT__ "PCView_MG"
3297233f9f0SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
3304b9ad928SBarry Smith {
3319dcbbd2bSBarry Smith   PC_MG          **mg = (PC_MG**)pc->data;
332dfbe8321SBarry Smith   PetscErrorCode ierr;
33379416396SBarry Smith   PetscInt       levels = mg[0]->levels,i;
33432077d6dSBarry Smith   PetscTruth     iascii;
3354b9ad928SBarry Smith 
3364b9ad928SBarry Smith   PetscFunctionBegin;
33732077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
33832077d6dSBarry Smith   if (iascii) {
339a8e6722dSMatthew 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);
340d792f543SMatthew Knepley     if (mg[0]->am == PC_MG_MULTIPLICATIVE) {
341d792f543SMatthew Knepley       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg[0]->cyclesperpcapply);CHKERRQ(ierr);
342d792f543SMatthew Knepley     }
343c2be2410SBarry Smith     if (mg[0]->galerkin) {
344c2be2410SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
345c2be2410SBarry Smith     }
3464b9ad928SBarry Smith     for (i=0; i<levels; i++) {
347b03c7568SBarry Smith       if (!i) {
348a8e6722dSMatthew 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);
349b03c7568SBarry Smith       } else {
350a8e6722dSMatthew Knepley         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg[i]->default_smoothd);CHKERRQ(ierr);
351b03c7568SBarry Smith       }
3524b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3534b9ad928SBarry Smith       ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr);
3544b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
355b03c7568SBarry Smith       if (i && mg[i]->smoothd == mg[i]->smoothu) {
3564b9ad928SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
357b03c7568SBarry Smith       } else if (i){
358a8e6722dSMatthew Knepley         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg[i]->default_smoothu);CHKERRQ(ierr);
3594b9ad928SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
3604b9ad928SBarry Smith         ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr);
3614b9ad928SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
3624b9ad928SBarry Smith       }
3634b9ad928SBarry Smith     }
3644b9ad928SBarry Smith   } else {
36579a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name);
3664b9ad928SBarry Smith   }
3674b9ad928SBarry Smith   PetscFunctionReturn(0);
3684b9ad928SBarry Smith }
3694b9ad928SBarry Smith 
3704b9ad928SBarry Smith /*
3714b9ad928SBarry Smith     Calls setup for the KSP on each level
3724b9ad928SBarry Smith */
3734b9ad928SBarry Smith #undef __FUNCT__
3744b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_MG"
3757233f9f0SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
3764b9ad928SBarry Smith {
3779dcbbd2bSBarry Smith   PC_MG                   **mg = (PC_MG**)pc->data;
378dfbe8321SBarry Smith   PetscErrorCode          ierr;
37979416396SBarry Smith   PetscInt                i,n = mg[0]->levels;
38077122347SBarry Smith   PC                      cpc,mpc;
381*90d69ab7SBarry Smith   PetscTruth              preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset;
38223d894e5SBarry Smith   PetscViewerASCIIMonitor ascii;
38323d894e5SBarry Smith   PetscViewer             viewer = PETSC_NULL;
3844b9ad928SBarry Smith   MPI_Comm                comm;
385c2be2410SBarry Smith   Mat                     dA,dB;
386c2be2410SBarry Smith   MatStructure            uflag;
3870a6bb862SBarry Smith   Vec                     tvec;
3884b9ad928SBarry Smith 
3894b9ad928SBarry Smith   PetscFunctionBegin;
390852665d3SBarry Smith 
39143fb2f97SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
39243fb2f97SBarry Smith   /* so use those from global PC */
39343fb2f97SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
394906ed7ccSBarry Smith   ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
39543fb2f97SBarry Smith   ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr);
39677122347SBarry Smith   ierr = KSPGetPC(mg[n-1]->smoothd,&mpc);CHKERRQ(ierr);
39777122347SBarry Smith   if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2))) {
398852665d3SBarry 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);
3991cfe3bddSBarry Smith     ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
400852665d3SBarry Smith   }
401852665d3SBarry Smith 
402852665d3SBarry Smith   if (mg[0]->galerkin) {
403852665d3SBarry Smith     Mat B;
404852665d3SBarry Smith     mg[0]->galerkinused = PETSC_TRUE;
405852665d3SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
406852665d3SBarry Smith     ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
407852665d3SBarry Smith     if (!pc->setupcalled) {
408852665d3SBarry Smith       for (i=n-2; i>-1; i--) {
409852665d3SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
410852665d3SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
411906ed7ccSBarry Smith 	if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
412852665d3SBarry Smith         dB   = B;
413852665d3SBarry Smith       }
414906ed7ccSBarry Smith       ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);
415852665d3SBarry Smith     } else {
416852665d3SBarry Smith       for (i=n-2; i>-1; i--) {
417906ed7ccSBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
418852665d3SBarry Smith         ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
419852665d3SBarry Smith         ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
420852665d3SBarry Smith         dB   = B;
421852665d3SBarry Smith       }
422852665d3SBarry Smith     }
423852665d3SBarry Smith   }
424852665d3SBarry Smith 
425958c9bccSBarry Smith   if (!pc->setupcalled) {
426*90d69ab7SBarry Smith     ierr = PetscOptionsGetTruth(0,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr);
4274b9ad928SBarry Smith 
428b03c7568SBarry Smith     for (i=0; i<n; i++) {
4294b9ad928SBarry Smith       if (monitor) {
4304b9ad928SBarry Smith         ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr);
43123d894e5SBarry Smith         ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
43223d894e5SBarry Smith         ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
4334b9ad928SBarry Smith       }
4344b9ad928SBarry Smith       ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr);
4354b9ad928SBarry Smith     }
436b03c7568SBarry Smith     for (i=1; i<n; i++) {
437a98fc643SBarry Smith       if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) {
4384b9ad928SBarry Smith         if (monitor) {
4394b9ad928SBarry Smith           ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr);
44023d894e5SBarry Smith           ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr);
44123d894e5SBarry Smith           ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
4424b9ad928SBarry Smith         }
4434b9ad928SBarry Smith         ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr);
4444b9ad928SBarry Smith       }
4454b9ad928SBarry Smith     }
446fccaa45eSBarry Smith     for (i=1; i<n; i++) {
4470cace4b0SBarry Smith       if (!mg[i]->residual) {
4480cace4b0SBarry Smith         Mat mat;
4490cace4b0SBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
4500cace4b0SBarry Smith         ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
4510cace4b0SBarry Smith       }
452fccaa45eSBarry Smith       if (mg[i]->restrct && !mg[i]->interpolate) {
453aea2a34eSBarry Smith         ierr = PCMGSetInterpolation(pc,i,mg[i]->restrct);CHKERRQ(ierr);
454fccaa45eSBarry Smith       }
455fccaa45eSBarry Smith       if (!mg[i]->restrct && mg[i]->interpolate) {
45697177400SBarry Smith         ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr);
457fccaa45eSBarry Smith       }
458fccaa45eSBarry Smith #if defined(PETSC_USE_DEBUG)
459fccaa45eSBarry Smith       if (!mg[i]->restrct || !mg[i]->interpolate) {
460fccaa45eSBarry Smith         SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i);
461fccaa45eSBarry Smith       }
462fccaa45eSBarry Smith #endif
463fccaa45eSBarry Smith     }
4640a6bb862SBarry Smith     for (i=0; i<n-1; i++) {
46537b0e6c0SBarry Smith       if (!mg[i]->b) {
466906ed7ccSBarry Smith         Vec *vec;
467906ed7ccSBarry Smith         ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
468906ed7ccSBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
4694fdc791cSSatish Balay         ierr = VecDestroy(*vec);CHKERRQ(ierr);
470906ed7ccSBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
47137b0e6c0SBarry Smith       }
4726ca4d86aSHong Zhang       if (!mg[i]->r && i) {
4730a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
47497177400SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
4750a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4760a6bb862SBarry Smith       }
4770a6bb862SBarry Smith       if (!mg[i]->x) {
4780a6bb862SBarry Smith         ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr);
47997177400SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
4800a6bb862SBarry Smith         ierr = VecDestroy(tvec);CHKERRQ(ierr);
4810a6bb862SBarry Smith       }
4820a6bb862SBarry Smith     }
483dfef70bdSHong Zhang     if (n != 1 && !mg[n-1]->r) {
484dfef70bdSHong Zhang       /* PCMGSetR() on the finest level if user did not supply it */
4850f77fa87SBarry Smith       Vec *vec;
4860f77fa87SBarry Smith       ierr = KSPGetVecs(mg[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
4870f77fa87SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
488797d2ccbSBarry Smith       ierr = VecDestroy(*vec);CHKERRQ(ierr);
4890f77fa87SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
4900f77fa87SBarry Smith     }
4914b9ad928SBarry Smith   }
4924b9ad928SBarry Smith 
493c2be2410SBarry Smith 
4944b9ad928SBarry Smith   for (i=1; i<n; i++) {
495b03c7568SBarry Smith     if (mg[i]->smoothu == mg[i]->smoothd) {
496b03c7568SBarry Smith       /* if doing only down then initial guess is zero */
4974b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
498b03c7568SBarry Smith     }
49932cf1786SBarry Smith     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5004b9ad928SBarry Smith     ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr);
50132cf1786SBarry Smith     if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5024b9ad928SBarry Smith   }
503b03c7568SBarry Smith   for (i=1; i<n; i++) {
5044b9ad928SBarry Smith     if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) {
505906ed7ccSBarry Smith       Mat          downmat,downpmat;
50697f1f81fSBarry Smith       MatStructure matflag;
507906ed7ccSBarry Smith       PetscTruth   opsset;
50897f1f81fSBarry Smith 
50997f1f81fSBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
510906ed7ccSBarry Smith       ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
511906ed7ccSBarry Smith       if (!opsset) {
512906ed7ccSBarry Smith         ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
51397f1f81fSBarry Smith         ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
51497f1f81fSBarry Smith       }
51597f1f81fSBarry Smith 
5164b9ad928SBarry Smith       ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
51732cf1786SBarry Smith       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5184b9ad928SBarry Smith       ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr);
51932cf1786SBarry Smith       if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5204b9ad928SBarry Smith     }
5214b9ad928SBarry Smith   }
5224b9ad928SBarry Smith 
5234b9ad928SBarry Smith   /*
5244b9ad928SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
5254b9ad928SBarry Smith   */
5264b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
5274b9ad928SBarry Smith   if (preonly) {
5284b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
5294b9ad928SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
53068eff7e6SBarry Smith     ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
53168eff7e6SBarry Smith     if (!lu && !redundant && !cholesky) {
5324b9ad928SBarry Smith       ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
5334b9ad928SBarry Smith     }
5344b9ad928SBarry Smith   }
5354b9ad928SBarry Smith 
536958c9bccSBarry Smith   if (!pc->setupcalled) {
5374b9ad928SBarry Smith     if (monitor) {
5384b9ad928SBarry Smith       ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr);
53923d894e5SBarry Smith       ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr);
54023d894e5SBarry Smith       ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr);
5414b9ad928SBarry Smith     }
5424b9ad928SBarry Smith     ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr);
5434b9ad928SBarry Smith   }
5444b9ad928SBarry Smith 
54532cf1786SBarry Smith   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5464b9ad928SBarry Smith   ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr);
54732cf1786SBarry Smith   if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
5484b9ad928SBarry Smith 
5494b9ad928SBarry Smith   /*
5506805f65bSBarry Smith      Dump the interpolation/restriction matrices plus the
5514b9ad928SBarry Smith    Jacobian/stiffness on each level. This allows Matlab users to
5526805f65bSBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
5536805f65bSBarry Smith 
5546805f65bSBarry Smith    Only support one or the other at the same time.
5556805f65bSBarry Smith   */
5566805f65bSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
557*90d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
5584b9ad928SBarry Smith   if (dump) {
5597adad957SLisandro Dalcin     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
5604b9ad928SBarry Smith   }
561*90d69ab7SBarry Smith   dump = PETSC_FALSE;
562c45a1595SBarry Smith #endif
563*90d69ab7SBarry Smith   ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
564c2be2410SBarry Smith   if (dump) {
5657adad957SLisandro Dalcin     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
5666805f65bSBarry Smith   }
5676805f65bSBarry Smith 
5686805f65bSBarry Smith   if (viewer) {
569c2be2410SBarry Smith     for (i=1; i<n; i++) {
5706805f65bSBarry Smith       ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr);
571c2be2410SBarry Smith     }
572c2be2410SBarry Smith     for (i=0; i<n; i++) {
573c2be2410SBarry Smith       ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr);
5746805f65bSBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
575c2be2410SBarry Smith     }
576c2be2410SBarry Smith   }
5774b9ad928SBarry Smith   PetscFunctionReturn(0);
5784b9ad928SBarry Smith }
5794b9ad928SBarry Smith 
5804b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
5814b9ad928SBarry Smith 
5824b9ad928SBarry Smith #undef __FUNCT__
5839dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
5844b9ad928SBarry Smith /*@C
58597177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
5864b9ad928SBarry Smith    Must be called before any other MG routine.
5874b9ad928SBarry Smith 
5884b9ad928SBarry Smith    Collective on PC
5894b9ad928SBarry Smith 
5904b9ad928SBarry Smith    Input Parameters:
5914b9ad928SBarry Smith +  pc - the preconditioner context
5924b9ad928SBarry Smith .  levels - the number of levels
5934b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
5944b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
5954b9ad928SBarry Smith 
5964b9ad928SBarry Smith    Level: intermediate
5974b9ad928SBarry Smith 
5984b9ad928SBarry Smith    Notes:
5994b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
6004b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
6014b9ad928SBarry Smith 
6024b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
6034b9ad928SBarry Smith 
60497177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
6054b9ad928SBarry Smith @*/
60697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
6074b9ad928SBarry Smith {
608dfbe8321SBarry Smith   PetscErrorCode ierr;
609ada7143aSSatish Balay   PC_MG          **mg=0;
6104b9ad928SBarry Smith 
6114b9ad928SBarry Smith   PetscFunctionBegin;
6124482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6134b9ad928SBarry Smith 
6144b9ad928SBarry Smith   if (pc->data) {
6151302d50aSBarry Smith     SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\
61697177400SBarry Smith     make sure that you call PCMGSetLevels() before KSPSetFromOptions()");
6174b9ad928SBarry Smith   }
6187adad957SLisandro Dalcin   ierr                     = PCMGCreate_Private(((PetscObject)pc)->comm,levels,pc,comms,&mg);CHKERRQ(ierr);
6199dcbbd2bSBarry Smith   mg[0]->am                = PC_MG_MULTIPLICATIVE;
6204b9ad928SBarry Smith   pc->data                 = (void*)mg;
6214b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
6224b9ad928SBarry Smith   PetscFunctionReturn(0);
6234b9ad928SBarry Smith }
6244b9ad928SBarry Smith 
6254b9ad928SBarry Smith #undef __FUNCT__
6269dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
6274b9ad928SBarry Smith /*@
62897177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
6294b9ad928SBarry Smith 
6304b9ad928SBarry Smith    Not Collective
6314b9ad928SBarry Smith 
6324b9ad928SBarry Smith    Input Parameter:
6334b9ad928SBarry Smith .  pc - the preconditioner context
6344b9ad928SBarry Smith 
6354b9ad928SBarry Smith    Output parameter:
6364b9ad928SBarry Smith .  levels - the number of levels
6374b9ad928SBarry Smith 
6384b9ad928SBarry Smith    Level: advanced
6394b9ad928SBarry Smith 
6404b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
6414b9ad928SBarry Smith 
64297177400SBarry Smith .seealso: PCMGSetLevels()
6434b9ad928SBarry Smith @*/
64497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels)
6454b9ad928SBarry Smith {
6469dcbbd2bSBarry Smith   PC_MG  **mg;
6474b9ad928SBarry Smith 
6484b9ad928SBarry Smith   PetscFunctionBegin;
6494482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6504482741eSBarry Smith   PetscValidIntPointer(levels,2);
6514b9ad928SBarry Smith 
6529dcbbd2bSBarry Smith   mg      = (PC_MG**)pc->data;
6534b9ad928SBarry Smith   *levels = mg[0]->levels;
6544b9ad928SBarry Smith   PetscFunctionReturn(0);
6554b9ad928SBarry Smith }
6564b9ad928SBarry Smith 
6574b9ad928SBarry Smith #undef __FUNCT__
6589dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
6594b9ad928SBarry Smith /*@
66097177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
6614b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
6624b9ad928SBarry Smith 
6634b9ad928SBarry Smith    Collective on PC
6644b9ad928SBarry Smith 
6654b9ad928SBarry Smith    Input Parameters:
6664b9ad928SBarry Smith +  pc - the preconditioner context
6679dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
6689dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
6694b9ad928SBarry Smith 
6704b9ad928SBarry Smith    Options Database Key:
6714b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
6724b9ad928SBarry Smith    additive, full, kaskade
6734b9ad928SBarry Smith 
6744b9ad928SBarry Smith    Level: advanced
6754b9ad928SBarry Smith 
6764b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
6774b9ad928SBarry Smith 
67897177400SBarry Smith .seealso: PCMGSetLevels()
6794b9ad928SBarry Smith @*/
6809dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form)
6814b9ad928SBarry Smith {
6829dcbbd2bSBarry Smith   PC_MG **mg;
6834b9ad928SBarry Smith 
6844b9ad928SBarry Smith   PetscFunctionBegin;
6854482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
6869dcbbd2bSBarry Smith   mg = (PC_MG**)pc->data;
6874b9ad928SBarry Smith 
6884b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
6894b9ad928SBarry Smith   mg[0]->am = form;
6909dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
6914b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
6924b9ad928SBarry Smith   PetscFunctionReturn(0);
6934b9ad928SBarry Smith }
6944b9ad928SBarry Smith 
6954b9ad928SBarry Smith #undef __FUNCT__
6960d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
6974b9ad928SBarry Smith /*@
6980d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
6994b9ad928SBarry Smith    complicated cycling.
7004b9ad928SBarry Smith 
7014b9ad928SBarry Smith    Collective on PC
7024b9ad928SBarry Smith 
7034b9ad928SBarry Smith    Input Parameters:
704c2be2410SBarry Smith +  pc - the multigrid context
7050d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
7064b9ad928SBarry Smith 
7074b9ad928SBarry Smith    Options Database Key:
7080d353602SBarry Smith $  -pc_mg_cycle_type v or w
7094b9ad928SBarry Smith 
7104b9ad928SBarry Smith    Level: advanced
7114b9ad928SBarry Smith 
7124b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7134b9ad928SBarry Smith 
7140d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
7154b9ad928SBarry Smith @*/
7160d353602SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n)
7174b9ad928SBarry Smith {
7189dcbbd2bSBarry Smith   PC_MG    **mg;
71979416396SBarry Smith   PetscInt i,levels;
7204b9ad928SBarry Smith 
7214b9ad928SBarry Smith   PetscFunctionBegin;
7224482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7239dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
7244b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7254b9ad928SBarry Smith   levels = mg[0]->levels;
7264b9ad928SBarry Smith 
7274b9ad928SBarry Smith   for (i=0; i<levels; i++) {
7284b9ad928SBarry Smith     mg[i]->cycles  = n;
7294b9ad928SBarry Smith   }
7304b9ad928SBarry Smith   PetscFunctionReturn(0);
7314b9ad928SBarry Smith }
7324b9ad928SBarry Smith 
7334b9ad928SBarry Smith #undef __FUNCT__
7348cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
7358cc2d5dfSBarry Smith /*@
7368cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
7378cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
7388cc2d5dfSBarry Smith 
7398cc2d5dfSBarry Smith    Collective on PC
7408cc2d5dfSBarry Smith 
7418cc2d5dfSBarry Smith    Input Parameters:
7428cc2d5dfSBarry Smith +  pc - the multigrid context
7438cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
7448cc2d5dfSBarry Smith 
7458cc2d5dfSBarry Smith    Options Database Key:
7468cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
7478cc2d5dfSBarry Smith 
7488cc2d5dfSBarry Smith    Level: advanced
7498cc2d5dfSBarry Smith 
7508cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
7518cc2d5dfSBarry Smith 
7528cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
7538cc2d5dfSBarry Smith 
7548cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
7558cc2d5dfSBarry Smith @*/
7568cc2d5dfSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
7578cc2d5dfSBarry Smith {
7588cc2d5dfSBarry Smith   PC_MG    **mg;
7598cc2d5dfSBarry Smith   PetscInt i,levels;
7608cc2d5dfSBarry Smith 
7618cc2d5dfSBarry Smith   PetscFunctionBegin;
7628cc2d5dfSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
7638cc2d5dfSBarry Smith   mg     = (PC_MG**)pc->data;
7648cc2d5dfSBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
7658cc2d5dfSBarry Smith   levels = mg[0]->levels;
7668cc2d5dfSBarry Smith 
7678cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
7688cc2d5dfSBarry Smith     mg[i]->cyclesperpcapply  = n;
7698cc2d5dfSBarry Smith   }
7708cc2d5dfSBarry Smith   PetscFunctionReturn(0);
7718cc2d5dfSBarry Smith }
7728cc2d5dfSBarry Smith 
7738cc2d5dfSBarry Smith #undef __FUNCT__
7749dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
775c2be2410SBarry Smith /*@
77697177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
777c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
778c2be2410SBarry Smith 
779c2be2410SBarry Smith    Collective on PC
780c2be2410SBarry Smith 
781c2be2410SBarry Smith    Input Parameters:
7823fc8bf9cSBarry Smith .  pc - the multigrid context
783c2be2410SBarry Smith 
784c2be2410SBarry Smith    Options Database Key:
785c2be2410SBarry Smith $  -pc_mg_galerkin
786c2be2410SBarry Smith 
787c2be2410SBarry Smith    Level: intermediate
788c2be2410SBarry Smith 
789c2be2410SBarry Smith .keywords: MG, set, Galerkin
790c2be2410SBarry Smith 
7913fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
7923fc8bf9cSBarry Smith 
793c2be2410SBarry Smith @*/
79497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc)
795c2be2410SBarry Smith {
7969dcbbd2bSBarry Smith   PC_MG    **mg;
797c2be2410SBarry Smith   PetscInt i,levels;
798c2be2410SBarry Smith 
799c2be2410SBarry Smith   PetscFunctionBegin;
800c2be2410SBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8019dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
802c2be2410SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
803c2be2410SBarry Smith   levels = mg[0]->levels;
804c2be2410SBarry Smith 
805c2be2410SBarry Smith   for (i=0; i<levels; i++) {
806c2be2410SBarry Smith     mg[i]->galerkin = PETSC_TRUE;
807c2be2410SBarry Smith   }
808c2be2410SBarry Smith   PetscFunctionReturn(0);
809c2be2410SBarry Smith }
810c2be2410SBarry Smith 
811c2be2410SBarry Smith #undef __FUNCT__
8123fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
8133fc8bf9cSBarry Smith /*@
8143fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
8153fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
8163fc8bf9cSBarry Smith 
8173fc8bf9cSBarry Smith    Not Collective
8183fc8bf9cSBarry Smith 
8193fc8bf9cSBarry Smith    Input Parameter:
8203fc8bf9cSBarry Smith .  pc - the multigrid context
8213fc8bf9cSBarry Smith 
8223fc8bf9cSBarry Smith    Output Parameter:
8233fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
8243fc8bf9cSBarry Smith 
8253fc8bf9cSBarry Smith    Options Database Key:
8263fc8bf9cSBarry Smith $  -pc_mg_galerkin
8273fc8bf9cSBarry Smith 
8283fc8bf9cSBarry Smith    Level: intermediate
8293fc8bf9cSBarry Smith 
8303fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
8313fc8bf9cSBarry Smith 
8323fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
8333fc8bf9cSBarry Smith 
8343fc8bf9cSBarry Smith @*/
8353fc8bf9cSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin)
8363fc8bf9cSBarry Smith {
8373fc8bf9cSBarry Smith   PC_MG    **mg;
8383fc8bf9cSBarry Smith 
8393fc8bf9cSBarry Smith   PetscFunctionBegin;
8403fc8bf9cSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8413fc8bf9cSBarry Smith   mg     = (PC_MG**)pc->data;
8423fc8bf9cSBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
8433fc8bf9cSBarry Smith   *galerkin = mg[0]->galerkin;
8443fc8bf9cSBarry Smith   PetscFunctionReturn(0);
8453fc8bf9cSBarry Smith }
8463fc8bf9cSBarry Smith 
8473fc8bf9cSBarry Smith #undef __FUNCT__
8489dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
8494b9ad928SBarry Smith /*@
85097177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
85197177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
8524b9ad928SBarry Smith    pre-smoothing steps on different levels.
8534b9ad928SBarry Smith 
8544b9ad928SBarry Smith    Collective on PC
8554b9ad928SBarry Smith 
8564b9ad928SBarry Smith    Input Parameters:
8574b9ad928SBarry Smith +  mg - the multigrid context
8584b9ad928SBarry Smith -  n - the number of smoothing steps
8594b9ad928SBarry Smith 
8604b9ad928SBarry Smith    Options Database Key:
8614b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
8624b9ad928SBarry Smith 
8634b9ad928SBarry Smith    Level: advanced
8644b9ad928SBarry Smith 
8654b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
8664b9ad928SBarry Smith 
86797177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
8684b9ad928SBarry Smith @*/
86997177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n)
8704b9ad928SBarry Smith {
8719dcbbd2bSBarry Smith   PC_MG          **mg;
8726849ba73SBarry Smith   PetscErrorCode ierr;
87379416396SBarry Smith   PetscInt       i,levels;
8744b9ad928SBarry Smith 
8754b9ad928SBarry Smith   PetscFunctionBegin;
8764482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
8779dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
8784b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
8794b9ad928SBarry Smith   levels = mg[0]->levels;
8804b9ad928SBarry Smith 
881b05257ddSBarry Smith   for (i=1; i<levels; i++) {
8824b9ad928SBarry Smith     /* make sure smoother up and down are different */
88397177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
8844b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
8854b9ad928SBarry Smith     mg[i]->default_smoothd = n;
8864b9ad928SBarry Smith   }
8874b9ad928SBarry Smith   PetscFunctionReturn(0);
8884b9ad928SBarry Smith }
8894b9ad928SBarry Smith 
8904b9ad928SBarry Smith #undef __FUNCT__
8919dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
8924b9ad928SBarry Smith /*@
89397177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
89497177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
8954b9ad928SBarry Smith    post-smoothing steps on different levels.
8964b9ad928SBarry Smith 
8974b9ad928SBarry Smith    Collective on PC
8984b9ad928SBarry Smith 
8994b9ad928SBarry Smith    Input Parameters:
9004b9ad928SBarry Smith +  mg - the multigrid context
9014b9ad928SBarry Smith -  n - the number of smoothing steps
9024b9ad928SBarry Smith 
9034b9ad928SBarry Smith    Options Database Key:
9044b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
9054b9ad928SBarry Smith 
9064b9ad928SBarry Smith    Level: advanced
9074b9ad928SBarry Smith 
9084b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
909a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
9104b9ad928SBarry Smith 
9114b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
9124b9ad928SBarry Smith 
91397177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
9144b9ad928SBarry Smith @*/
91597177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n)
9164b9ad928SBarry Smith {
9179dcbbd2bSBarry Smith   PC_MG          **mg;
9186849ba73SBarry Smith   PetscErrorCode ierr;
91979416396SBarry Smith   PetscInt       i,levels;
9204b9ad928SBarry Smith 
9214b9ad928SBarry Smith   PetscFunctionBegin;
9224482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
9239dcbbd2bSBarry Smith   mg     = (PC_MG**)pc->data;
9244b9ad928SBarry Smith   if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
9254b9ad928SBarry Smith   levels = mg[0]->levels;
9264b9ad928SBarry Smith 
9274b9ad928SBarry Smith   for (i=1; i<levels; i++) {
9284b9ad928SBarry Smith     /* make sure smoother up and down are different */
92997177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
9304b9ad928SBarry Smith     ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
9314b9ad928SBarry Smith     mg[i]->default_smoothu = n;
9324b9ad928SBarry Smith   }
9334b9ad928SBarry Smith   PetscFunctionReturn(0);
9344b9ad928SBarry Smith }
9354b9ad928SBarry Smith 
9364b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
9374b9ad928SBarry Smith 
9383b09bd56SBarry Smith /*MC
939ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
9403b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
9413b09bd56SBarry Smith 
9423b09bd56SBarry Smith    Options Database Keys:
9433b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
9440d353602SBarry Smith .  -pc_mg_cycles v or w
94579416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
9463b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
9473b09bd56SBarry Smith .  -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default
9483b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
9493b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
95068eff7e6SBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators
9513b09bd56SBarry Smith -  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
9523b09bd56SBarry Smith                         to the Socket viewer for reading from Matlab.
9533b09bd56SBarry Smith 
9543b09bd56SBarry Smith    Notes:
9553b09bd56SBarry Smith 
9563b09bd56SBarry Smith    Level: intermediate
9573b09bd56SBarry Smith 
9588f87f92bSBarry Smith    Concepts: multigrid/multilevel
9593b09bd56SBarry Smith 
9603b09bd56SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
9610d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
96297177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
96397177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
9640d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
9653b09bd56SBarry Smith M*/
9663b09bd56SBarry Smith 
9674b9ad928SBarry Smith EXTERN_C_BEGIN
9684b9ad928SBarry Smith #undef __FUNCT__
9694b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
970dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc)
9714b9ad928SBarry Smith {
9724b9ad928SBarry Smith   PetscFunctionBegin;
9734b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
9744b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
9754b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
9764b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
9774b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
9784b9ad928SBarry Smith 
9794b9ad928SBarry Smith   pc->data                = (void*)0;
9804b9ad928SBarry Smith   PetscFunctionReturn(0);
9814b9ad928SBarry Smith }
9824b9ad928SBarry Smith EXTERN_C_END
983