xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 9dbfc18722704d6a1497b9a29671ad8f32ee6572)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h>                    /*I "petscpcmg.h" I*/
64b9ad928SBarry Smith 
74b9ad928SBarry Smith 
84b9ad928SBarry Smith #undef __FUNCT__
99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private"
1031567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
114b9ad928SBarry Smith {
1231567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1331567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
146849ba73SBarry Smith   PetscErrorCode ierr;
1531567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
1863e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
1931567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
2063e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2131567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2263e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2331567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2463e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
254b9ad928SBarry Smith 
264b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
2731567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
284b9ad928SBarry Smith       PetscReal rnorm;
2931567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
304b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3170441072SBarry Smith         if (rnorm < mg->abstol) {
324d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
331e2582c4SBarry Smith           ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr);
344b9ad928SBarry Smith         } else {
354d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
361e2582c4SBarry 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);
374b9ad928SBarry Smith         }
384b9ad928SBarry Smith         PetscFunctionReturn(0);
394b9ad928SBarry Smith       }
404b9ad928SBarry Smith     }
414b9ad928SBarry Smith 
4231567311SBarry Smith     mgc = *(mglevelsin - 1);
4363e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4431567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
4563e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
46efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
474b9ad928SBarry Smith     while (cycles--) {
4831567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
494b9ad928SBarry Smith     }
5063e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5131567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5263e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5363e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5431567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
5563e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
564b9ad928SBarry Smith   }
574b9ad928SBarry Smith   PetscFunctionReturn(0);
584b9ad928SBarry Smith }
594b9ad928SBarry Smith 
604b9ad928SBarry Smith #undef __FUNCT__
614b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG"
62ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool  zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason)
634b9ad928SBarry Smith {
64f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
65f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
66dfbe8321SBarry Smith   PetscErrorCode ierr;
67f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
684b9ad928SBarry Smith 
694b9ad928SBarry Smith   PetscFunctionBegin;
70f3fbd535SBarry Smith   mglevels[levels-1]->b    = b;
71f3fbd535SBarry Smith   mglevels[levels-1]->x    = x;
724b9ad928SBarry Smith 
7331567311SBarry Smith   mg->rtol = rtol;
7431567311SBarry Smith   mg->abstol = abstol;
7531567311SBarry Smith   mg->dtol = dtol;
764b9ad928SBarry Smith   if (rtol) {
774b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
784b9ad928SBarry Smith     PetscReal rnorm;
797319c654SBarry Smith     if (zeroguess) {
807319c654SBarry Smith       ierr               = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
817319c654SBarry Smith     } else {
82f3fbd535SBarry Smith       ierr               = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
834b9ad928SBarry Smith       ierr               = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
847319c654SBarry Smith     }
8531567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
8670441072SBarry Smith   } else if (abstol) {
8731567311SBarry Smith     mg->ttol = abstol;
884b9ad928SBarry Smith   } else {
8931567311SBarry Smith     mg->ttol = 0.0;
904b9ad928SBarry Smith   }
914b9ad928SBarry Smith 
924d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
93336babb1SJed Brown      stop prematurely due to small residual */
944d0a8057SBarry Smith   for (i=1; i<levels; i++) {
95f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
96f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
97f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
984b9ad928SBarry Smith     }
994d0a8057SBarry Smith   }
1004d0a8057SBarry Smith 
1014d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1024d0a8057SBarry Smith   for (i=0; i<its; i++) {
103f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1044d0a8057SBarry Smith     if (*reason) break;
1054d0a8057SBarry Smith   }
1064d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1074d0a8057SBarry Smith   *outits = i;
1084b9ad928SBarry Smith   PetscFunctionReturn(0);
1094b9ad928SBarry Smith }
1104b9ad928SBarry Smith 
1114b9ad928SBarry Smith #undef __FUNCT__
1123aeaf226SBarry Smith #define __FUNCT__ "PCReset_MG"
1133aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1143aeaf226SBarry Smith {
1153aeaf226SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1163aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1173aeaf226SBarry Smith   PetscErrorCode ierr;
1183aeaf226SBarry Smith   PetscInt       i,n;
1193aeaf226SBarry Smith 
1203aeaf226SBarry Smith   PetscFunctionBegin;
1213aeaf226SBarry Smith   if (mglevels) {
1223aeaf226SBarry Smith     n = mglevels[0]->levels;
1233aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1243aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1253aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1263aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
1273aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1283aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
12973250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1303aeaf226SBarry Smith     }
1313aeaf226SBarry Smith 
1323aeaf226SBarry Smith     for (i=0; i<n; i++) {
1333aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1343aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1353aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1363aeaf226SBarry Smith       }
1373aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1383aeaf226SBarry Smith     }
1393aeaf226SBarry Smith   }
1403aeaf226SBarry Smith   PetscFunctionReturn(0);
1413aeaf226SBarry Smith }
1423aeaf226SBarry Smith 
1433aeaf226SBarry Smith #undef __FUNCT__
1449dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels"
1454b9ad928SBarry Smith /*@C
14697177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1474b9ad928SBarry Smith    Must be called before any other MG routine.
1484b9ad928SBarry Smith 
149ad4df100SBarry Smith    Logically Collective on PC
1504b9ad928SBarry Smith 
1514b9ad928SBarry Smith    Input Parameters:
1524b9ad928SBarry Smith +  pc - the preconditioner context
1534b9ad928SBarry Smith .  levels - the number of levels
1544b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1554b9ad928SBarry Smith            on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith    Level: intermediate
1584b9ad928SBarry Smith 
1594b9ad928SBarry Smith    Notes:
1604b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1614b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1624b9ad928SBarry Smith 
1634b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1644b9ad928SBarry Smith 
16597177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1664b9ad928SBarry Smith @*/
1677087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1684b9ad928SBarry Smith {
169dfbe8321SBarry Smith   PetscErrorCode ierr;
170f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
171f3fbd535SBarry Smith   MPI_Comm       comm = ((PetscObject)pc)->comm;
1723aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
173f3fbd535SBarry Smith   PetscInt       i;
174f3fbd535SBarry Smith   PetscMPIInt    size;
175f3fbd535SBarry Smith   const char     *prefix;
176f3fbd535SBarry Smith   PC             ipc;
1773aeaf226SBarry Smith   PetscInt       n;
1784b9ad928SBarry Smith 
1794b9ad928SBarry Smith   PetscFunctionBegin;
1800700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
181548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
182548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1833aeaf226SBarry Smith   if (mglevels) {
1843aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
1853aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
1863aeaf226SBarry Smith     n = mglevels[0]->levels;
1873aeaf226SBarry Smith     for (i=0; i<n; i++) {
1883aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1893aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
1903aeaf226SBarry Smith       }
1913aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
1923aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
1933aeaf226SBarry Smith     }
1943aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
1953aeaf226SBarry Smith   }
196f3fbd535SBarry Smith 
197f3fbd535SBarry Smith   mg->nlevels      = levels;
198f3fbd535SBarry Smith 
199f3fbd535SBarry Smith   ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr);
200f3fbd535SBarry Smith   ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
201f3fbd535SBarry Smith 
202f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
203f3fbd535SBarry Smith 
204a9db87a2SMatthew G Knepley   mg->stageApply = 0;
205f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
206f3fbd535SBarry Smith     ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr);
207f3fbd535SBarry Smith     mglevels[i]->level           = i;
208f3fbd535SBarry Smith     mglevels[i]->levels          = levels;
209f3fbd535SBarry Smith     mglevels[i]->cycles          = PC_MG_CYCLE_V;
210336babb1SJed Brown     mg->default_smoothu = 2;
211336babb1SJed Brown     mg->default_smoothd = 2;
21263e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
21363e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
21463e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
21563e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
216f3fbd535SBarry Smith 
217f3fbd535SBarry Smith     if (comms) comm = comms[i];
218f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
219336babb1SJed Brown     ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
220336babb1SJed Brown     ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
221336babb1SJed Brown     ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
222336babb1SJed Brown     ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
223336babb1SJed Brown     ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
224f3fbd535SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
225336babb1SJed Brown     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i?mg->default_smoothd:1);CHKERRQ(ierr);
226f3fbd535SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
227f3fbd535SBarry Smith 
228f3fbd535SBarry Smith     /* do special stuff for coarse grid */
229f3fbd535SBarry Smith     if (!i && levels > 1) {
230f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
231f3fbd535SBarry Smith 
232*9dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
233f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
234f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
235f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
236f3fbd535SBarry Smith       if (size > 1) {
237*9dbfc187SHong Zhang         KSP innerksp;
238*9dbfc187SHong Zhang         PC  innerpc;
239f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
240*9dbfc187SHong Zhang         ierr = PCRedundantGetKSP(ipc,&innerksp);CHKERRQ(ierr);
241*9dbfc187SHong Zhang         ierr = KSPGetPC(innerksp,&innerpc);CHKERRQ(ierr);
242*9dbfc187SHong Zhang         ierr = PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
243f3fbd535SBarry Smith       } else {
244f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
245*9dbfc187SHong Zhang         ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
246f3fbd535SBarry Smith       }
247f3fbd535SBarry Smith     } else {
248f3fbd535SBarry Smith       char tprefix[128];
249f3fbd535SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
250f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
251f3fbd535SBarry Smith     }
252f3fbd535SBarry Smith     ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr);
253f3fbd535SBarry Smith     mglevels[i]->smoothu    = mglevels[i]->smoothd;
25431567311SBarry Smith     mg->rtol                = 0.0;
25531567311SBarry Smith     mg->abstol              = 0.0;
25631567311SBarry Smith     mg->dtol                = 0.0;
25731567311SBarry Smith     mg->ttol                = 0.0;
25831567311SBarry Smith     mg->cyclesperpcapply    = 1;
259f3fbd535SBarry Smith   }
26031567311SBarry Smith   mg->am          = PC_MG_MULTIPLICATIVE;
261f3fbd535SBarry Smith   mg->levels      = mglevels;
2624b9ad928SBarry Smith   pc->ops->applyrichardson = PCApplyRichardson_MG;
2634b9ad928SBarry Smith   PetscFunctionReturn(0);
2644b9ad928SBarry Smith }
2654b9ad928SBarry Smith 
266c07bf074SBarry Smith 
267c07bf074SBarry Smith #undef __FUNCT__
268c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG"
269c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
270c07bf074SBarry Smith {
271c07bf074SBarry Smith   PetscErrorCode ierr;
272a06653b4SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
273a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
274a06653b4SBarry Smith   PetscInt       i,n;
275c07bf074SBarry Smith 
276c07bf074SBarry Smith   PetscFunctionBegin;
277a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
278a06653b4SBarry Smith   if (mglevels) {
279a06653b4SBarry Smith     n = mglevels[0]->levels;
280a06653b4SBarry Smith     for (i=0; i<n; i++) {
281a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2826bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
283a06653b4SBarry Smith       }
2846bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
285a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
286a06653b4SBarry Smith     }
287a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
288a06653b4SBarry Smith   }
289c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
290f3fbd535SBarry Smith   PetscFunctionReturn(0);
291f3fbd535SBarry Smith }
292f3fbd535SBarry Smith 
293f3fbd535SBarry Smith 
294f3fbd535SBarry Smith 
29509573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
29609573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
29709573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
298f3fbd535SBarry Smith 
299f3fbd535SBarry Smith /*
300f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
301f3fbd535SBarry Smith              or full cycle of multigrid.
302f3fbd535SBarry Smith 
303f3fbd535SBarry Smith   Note:
304f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
305f3fbd535SBarry Smith */
306f3fbd535SBarry Smith #undef __FUNCT__
307f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG"
308f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
309f3fbd535SBarry Smith {
310f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
311f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
312f3fbd535SBarry Smith   PetscErrorCode ierr;
313f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
314f3fbd535SBarry Smith 
315f3fbd535SBarry Smith   PetscFunctionBegin;
316a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
317e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
318298cc208SBarry Smith   for (i=0; i<levels; i++) {
319e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
320e1d8e5deSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
321298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
322e1d8e5deSBarry Smith     }
323e1d8e5deSBarry Smith   }
324e1d8e5deSBarry Smith 
325f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
326f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
32731567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
328f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
32931567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
330f3fbd535SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr);
331f3fbd535SBarry Smith     }
332f3fbd535SBarry Smith   }
33331567311SBarry Smith   else if (mg->am == PC_MG_ADDITIVE) {
33431567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
335f3fbd535SBarry Smith   }
33631567311SBarry Smith   else if (mg->am == PC_MG_KASKADE) {
33731567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
338f3fbd535SBarry Smith   }
339f3fbd535SBarry Smith   else {
340f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
341f3fbd535SBarry Smith   }
342a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
343f3fbd535SBarry Smith   PetscFunctionReturn(0);
344f3fbd535SBarry Smith }
345f3fbd535SBarry Smith 
346f3fbd535SBarry Smith 
347f3fbd535SBarry Smith #undef __FUNCT__
348f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG"
349f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc)
350f3fbd535SBarry Smith {
351f3fbd535SBarry Smith   PetscErrorCode ierr;
352f3fbd535SBarry Smith   PetscInt       m,levels = 1,cycles;
353c91913d3SJed Brown   PetscBool      flg,set;
354f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
355f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
356f3fbd535SBarry Smith   PCMGType       mgtype;
357f3fbd535SBarry Smith   PCMGCycleType  mgctype;
358f3fbd535SBarry Smith 
359f3fbd535SBarry Smith   PetscFunctionBegin;
360f3fbd535SBarry Smith   ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr);
3613aeaf226SBarry Smith     if (!mg->levels) {
362f3fbd535SBarry Smith       ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
363eb3f98d2SBarry Smith       if (!flg && pc->dm) {
364eb3f98d2SBarry Smith         ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
365eb3f98d2SBarry Smith         levels++;
3663aeaf226SBarry Smith         mg->usedmfornumberoflevels = PETSC_TRUE;
367eb3f98d2SBarry Smith       }
368f3fbd535SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
369f3fbd535SBarry Smith     }
3703aeaf226SBarry Smith     mglevels = mg->levels;
3713aeaf226SBarry Smith 
372f3fbd535SBarry Smith     mgctype = (PCMGCycleType) mglevels[0]->cycles;
373f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
374f3fbd535SBarry Smith     if (flg) {
375f3fbd535SBarry Smith       ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
376f3fbd535SBarry Smith     };
377f3fbd535SBarry Smith     flg  = PETSC_FALSE;
378c91913d3SJed Brown     ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);CHKERRQ(ierr);
379c91913d3SJed Brown     if (set) {
380c91913d3SJed Brown       ierr = PCMGSetGalerkin(pc,flg);CHKERRQ(ierr);
381f3fbd535SBarry Smith     }
382f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
383f3fbd535SBarry Smith     if (flg) {
384f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
385f3fbd535SBarry Smith     }
386f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
387f3fbd535SBarry Smith     if (flg) {
388f3fbd535SBarry Smith       ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
389f3fbd535SBarry Smith     }
39031567311SBarry Smith     mgtype = mg->am;
391f3fbd535SBarry Smith     ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
392f3fbd535SBarry Smith     if (flg) {
393f3fbd535SBarry Smith       ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
394f3fbd535SBarry Smith     }
39531567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
39631567311SBarry Smith       ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
397f3fbd535SBarry Smith       if (flg) {
398f3fbd535SBarry Smith         ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
399f3fbd535SBarry Smith       }
400f3fbd535SBarry Smith     }
401f3fbd535SBarry Smith     flg  = PETSC_FALSE;
402acfcf0e5SJed Brown     ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
403f3fbd535SBarry Smith     if (flg) {
404f3fbd535SBarry Smith       PetscInt i;
405f3fbd535SBarry Smith       char     eventname[128];
40663e6d426SJed Brown       if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
407f3fbd535SBarry Smith       levels = mglevels[0]->levels;
408f3fbd535SBarry Smith       for (i=0; i<levels; i++) {
409f3fbd535SBarry Smith         sprintf(eventname,"MGSetup Level %d",(int)i);
41063e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
411f3fbd535SBarry Smith         sprintf(eventname,"MGSmooth Level %d",(int)i);
41263e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
413f3fbd535SBarry Smith         if (i) {
414f3fbd535SBarry Smith           sprintf(eventname,"MGResid Level %d",(int)i);
41563e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
416f3fbd535SBarry Smith           sprintf(eventname,"MGInterp Level %d",(int)i);
41763e6d426SJed Brown           ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
418f3fbd535SBarry Smith         }
419f3fbd535SBarry Smith       }
420eec35531SMatthew G Knepley 
421ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
422eec35531SMatthew G Knepley       {
423eec35531SMatthew G Knepley         const char   *sname = "MG Apply";
424eec35531SMatthew G Knepley         PetscStageLog stageLog;
425eec35531SMatthew G Knepley         PetscInt      st;
426eec35531SMatthew G Knepley 
427eec35531SMatthew G Knepley         PetscFunctionBegin;
428eec35531SMatthew G Knepley         ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
429eec35531SMatthew G Knepley         for (st = 0; st < stageLog->numStages; ++st) {
430eec35531SMatthew G Knepley           PetscBool same;
431eec35531SMatthew G Knepley 
432eec35531SMatthew G Knepley           ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
433eec35531SMatthew G Knepley           if (same) {mg->stageApply = st;}
434eec35531SMatthew G Knepley         }
435eec35531SMatthew G Knepley         if (!mg->stageApply) {
436eec35531SMatthew G Knepley           ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
437eec35531SMatthew G Knepley         }
438eec35531SMatthew G Knepley       }
439ec60ca73SSatish Balay #endif
440f3fbd535SBarry Smith     }
441f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
442f3fbd535SBarry Smith   PetscFunctionReturn(0);
443f3fbd535SBarry Smith }
444f3fbd535SBarry Smith 
4456a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4466a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
447f3fbd535SBarry Smith 
448f3fbd535SBarry Smith #undef __FUNCT__
449f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG"
450f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
451f3fbd535SBarry Smith {
452f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
453f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
454f3fbd535SBarry Smith   PetscErrorCode ierr;
455f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
456219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
457f3fbd535SBarry Smith 
458f3fbd535SBarry Smith   PetscFunctionBegin;
459251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4605b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
461219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
462f3fbd535SBarry Smith   if (iascii) {
46331567311SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,(mglevels[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr);
46431567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
46531567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
466f3fbd535SBarry Smith     }
467218a07d4SBarry Smith     if (mg->galerkin) {
468f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4694f66f45eSBarry Smith     } else {
4704f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
471f3fbd535SBarry Smith     }
472f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
473f3fbd535SBarry Smith       if (!i) {
47489cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
475f3fbd535SBarry Smith       } else {
47689cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
477f3fbd535SBarry Smith       }
478f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
479f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
480f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
481f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
482f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
483f3fbd535SBarry Smith       } else if (i) {
48489cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
485f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
486f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
487f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
488f3fbd535SBarry Smith       }
489f3fbd535SBarry Smith     }
4905b0b0462SBarry Smith   } else if (isbinary) {
4915b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
4925b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
4935b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
4945b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
4955b0b0462SBarry Smith       }
4965b0b0462SBarry Smith     }
497219b1045SBarry Smith   } else if (isdraw) {
498219b1045SBarry Smith     PetscDraw draw;
499b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
500219b1045SBarry Smith     ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
501219b1045SBarry Smith     ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
502219b1045SBarry Smith     ierr = PetscDrawStringGetSize(draw,PETSC_NULL,&th);CHKERRQ(ierr);
503219b1045SBarry Smith     bottom = y - th;
504219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
505b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
506219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
507219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
508219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
509b4375e8dSPeter Brune       } else {
510b4375e8dSPeter Brune         w = 0.5*PetscMin(1.0-x,x);
511b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
512b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
513b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
514b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
515b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
516b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
517b4375e8dSPeter Brune       }
5181cd381d6SBarry Smith       ierr = PetscDrawGetBoundingBox(draw,PETSC_NULL,&bottom,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
5191cd381d6SBarry Smith       bottom -= th;
520219b1045SBarry Smith     }
5215b0b0462SBarry Smith   }
522f3fbd535SBarry Smith   PetscFunctionReturn(0);
523f3fbd535SBarry Smith }
524f3fbd535SBarry Smith 
525b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>
526b45d2f2cSJed Brown #include <petsc-private/kspimpl.h>
527cab2e9ccSBarry Smith 
528f3fbd535SBarry Smith /*
529f3fbd535SBarry Smith     Calls setup for the KSP on each level
530f3fbd535SBarry Smith */
531f3fbd535SBarry Smith #undef __FUNCT__
532f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG"
533f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
534f3fbd535SBarry Smith {
535f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
536f3fbd535SBarry Smith   PC_MG_Levels            **mglevels = mg->levels;
537f3fbd535SBarry Smith   PetscErrorCode          ierr;
538f3fbd535SBarry Smith   PetscInt                i,n = mglevels[0]->levels;
53998e04984SBarry Smith   PC                      cpc;
540649052a6SBarry Smith   PetscBool               preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset;
541f3fbd535SBarry Smith   Mat                     dA,dB;
542f3fbd535SBarry Smith   MatStructure            uflag;
543f3fbd535SBarry Smith   Vec                     tvec;
544218a07d4SBarry Smith   DM                      *dms;
545649052a6SBarry Smith   PetscViewer             viewer = 0;
546f3fbd535SBarry Smith 
547f3fbd535SBarry Smith   PetscFunctionBegin;
54801bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5493aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5503aeaf226SBarry Smith     PetscInt levels;
5513aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5523aeaf226SBarry Smith     levels++;
5533aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5543aeaf226SBarry Smith       ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr);
5553aeaf226SBarry Smith       n    = levels;
5563aeaf226SBarry Smith       ierr = PCSetFromOptions(pc);CHKERRQ(ierr);  /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5573aeaf226SBarry Smith       mglevels =  mg->levels;
5583aeaf226SBarry Smith     }
5593aeaf226SBarry Smith   }
56098e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5613aeaf226SBarry Smith 
562f3fbd535SBarry Smith 
563f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
564f3fbd535SBarry Smith   /* so use those from global PC */
565f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
566f3fbd535SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr);
56798e04984SBarry Smith   if (opsset) {
56898e04984SBarry Smith     Mat mmat;
56998e04984SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,PETSC_NULL,&mmat,PETSC_NULL);CHKERRQ(ierr);
57098e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
57198e04984SBarry Smith   }
57298e04984SBarry Smith   if (!opsset) {
573f3fbd535SBarry 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);
574f3fbd535SBarry Smith     ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
575f3fbd535SBarry Smith   }
576f3fbd535SBarry Smith 
57701bc414fSDmitry Karpeev   /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */
578ce4cda84SJed Brown   if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
5792d2b81a6SBarry Smith     /* construct the interpolation from the DMs */
5802e499ae9SBarry Smith     Mat p;
58173250ac0SBarry Smith     Vec rscale;
582218a07d4SBarry Smith     ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr);
583218a07d4SBarry Smith     dms[n-1] = pc->dm;
584218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
585942e3340SBarry Smith       DMKSP kdm;
5862ee06e3bSJed Brown       ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);
5873c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
5883c0c59f3SBarry Smith       if (mg->galerkin) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
589942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
590d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
591d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
5925358d0d4SBarry Smith       kdm->ops->computerhs = PETSC_NULL;
593fe86f630SJed Brown       kdm->rhsctx          = PETSC_NULL;
59424c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
595e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
5962d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
59700ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
59873250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
5996bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
600218a07d4SBarry Smith       }
60124c3aa18SBarry Smith     }
6022d2b81a6SBarry Smith 
603218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
6046bf464f9SBarry Smith       ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);
605218a07d4SBarry Smith     }
6062d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
607ce4cda84SJed Brown   }
608cab2e9ccSBarry Smith 
609ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
610ce4cda84SJed Brown     /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
611cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
612cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
613218a07d4SBarry Smith   }
614218a07d4SBarry Smith 
615c91913d3SJed Brown   if (mg->galerkin == 1) {
616f3fbd535SBarry Smith     Mat B;
617f3fbd535SBarry Smith     /* currently only handle case where mat and pmat are the same on coarser levels */
618f3fbd535SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr);
619f3fbd535SBarry Smith     if (!pc->setupcalled) {
620f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
621f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr);
622f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
623f3fbd535SBarry Smith         if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
624f3fbd535SBarry Smith         dB   = B;
625f3fbd535SBarry Smith       }
626cd9507b2SBarry Smith       if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);}
627f3fbd535SBarry Smith     } else {
628f3fbd535SBarry Smith       for (i=n-2; i>-1; i--) {
629f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr);
630f3fbd535SBarry Smith         ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr);
631f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr);
632f3fbd535SBarry Smith         dB   = B;
633f3fbd535SBarry Smith       }
634f3fbd535SBarry Smith     }
635ce4cda84SJed Brown   } else if (!mg->galerkin && pc->dm && pc->dm->x) {
636cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
637cab2e9ccSBarry Smith     for (i=n-2;i>-1; i--) {
638c88c5224SJed Brown       Mat R;
639c88c5224SJed Brown       Vec rscale;
640cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
641cab2e9ccSBarry Smith         Vec *vecs;
642cab2e9ccSBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,PETSC_NULL);CHKERRQ(ierr);
643cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
644cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
645cab2e9ccSBarry Smith       }
646c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
647c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
648c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
649c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
650cab2e9ccSBarry Smith     }
651f3fbd535SBarry Smith   }
652ccceb50cSJed Brown   if (!mg->galerkin && pc->dm) {
653caa4e7f2SJed Brown     for (i=n-2;i>=0; i--) {
654ccceb50cSJed Brown       DM dmfine,dmcoarse;
655caa4e7f2SJed Brown       Mat Restrict,Inject;
656caa4e7f2SJed Brown       Vec rscale;
657ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
658ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
659caa4e7f2SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
660caa4e7f2SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
661caa4e7f2SJed Brown       Inject = PETSC_NULL;      /* Callback should create it if it needs Injection */
662caa4e7f2SJed Brown       ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
663caa4e7f2SJed Brown     }
664caa4e7f2SJed Brown   }
665f3fbd535SBarry Smith 
666f3fbd535SBarry Smith   if (!pc->setupcalled) {
667f3fbd535SBarry Smith     for (i=0; i<n; i++) {
668f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
669f3fbd535SBarry Smith     }
670f3fbd535SBarry Smith     for (i=1; i<n; i++) {
671f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
672f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
673f3fbd535SBarry Smith       }
674f3fbd535SBarry Smith     }
675f3fbd535SBarry Smith     for (i=1; i<n; i++) {
676c88c5224SJed Brown       ierr = PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);CHKERRQ(ierr);
677c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i,&mglevels[i]->restrct);CHKERRQ(ierr);
678f3fbd535SBarry Smith     }
679f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
680f3fbd535SBarry Smith       if (!mglevels[i]->b) {
681f3fbd535SBarry Smith         Vec *vec;
682f3fbd535SBarry Smith         ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
683f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
6846bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
685f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
686f3fbd535SBarry Smith       }
687f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
688f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
689f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
6906bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
691f3fbd535SBarry Smith       }
692f3fbd535SBarry Smith       if (!mglevels[i]->x) {
693f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
694f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
6956bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
696f3fbd535SBarry Smith       }
697f3fbd535SBarry Smith     }
698f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
699f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
700f3fbd535SBarry Smith       Vec *vec;
701f3fbd535SBarry Smith       ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr);
702f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
7036bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
704f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
705f3fbd535SBarry Smith     }
706f3fbd535SBarry Smith   }
707f3fbd535SBarry Smith 
70898e04984SBarry Smith   if (pc->dm) {
70998e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
71098e04984SBarry Smith     for (i=0; i<n-1; i++) {
71198e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
71298e04984SBarry Smith     }
71398e04984SBarry Smith   }
714f3fbd535SBarry Smith 
715f3fbd535SBarry Smith   for (i=1; i<n; i++) {
716f3fbd535SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
717f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
718f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
719f3fbd535SBarry Smith     }
72063e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
721f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
72263e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
723d42688cbSBarry Smith     if (!mglevels[i]->residual) {
724d42688cbSBarry Smith       Mat mat;
725d42688cbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr);
726d42688cbSBarry Smith       ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr);
727d42688cbSBarry Smith     }
728f3fbd535SBarry Smith   }
729f3fbd535SBarry Smith   for (i=1; i<n; i++) {
730f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
731f3fbd535SBarry Smith       Mat          downmat,downpmat;
732f3fbd535SBarry Smith       MatStructure matflag;
733f3fbd535SBarry Smith 
734f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
735f3fbd535SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr);
736f3fbd535SBarry Smith       if (!opsset) {
737f3fbd535SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr);
738f3fbd535SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr);
739f3fbd535SBarry Smith       }
740f3fbd535SBarry Smith 
741f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
74263e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
743f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
74463e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
745f3fbd535SBarry Smith     }
746f3fbd535SBarry Smith   }
747f3fbd535SBarry Smith 
748f3fbd535SBarry Smith   /*
749f3fbd535SBarry Smith       If coarse solver is not direct method then DO NOT USE preonly
750f3fbd535SBarry Smith   */
751251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr);
752f3fbd535SBarry Smith   if (preonly) {
753251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr);
754251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr);
755251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr);
756251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);CHKERRQ(ierr);
757fd1303e9SJungho Lee     if (!lu && !redundant && !cholesky && !svd) {
758f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr);
759f3fbd535SBarry Smith     }
760f3fbd535SBarry Smith   }
761f3fbd535SBarry Smith 
762f3fbd535SBarry Smith   if (!pc->setupcalled) {
763f3fbd535SBarry Smith     ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr);
764f3fbd535SBarry Smith   }
765f3fbd535SBarry Smith 
76663e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
767f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
76863e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
769f3fbd535SBarry Smith 
770f3fbd535SBarry Smith   /*
771f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
772e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
773f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
774f3fbd535SBarry Smith 
775f3fbd535SBarry Smith    Only support one or the other at the same time.
776f3fbd535SBarry Smith   */
777f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
778acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr);
779f3fbd535SBarry Smith   if (dump) {
780f3fbd535SBarry Smith     viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm);
781f3fbd535SBarry Smith   }
782f3fbd535SBarry Smith   dump = PETSC_FALSE;
783f3fbd535SBarry Smith #endif
784acfcf0e5SJed Brown   ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr);
785f3fbd535SBarry Smith   if (dump) {
78632c0f0efSBarry Smith 
787f3fbd535SBarry Smith     viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm);
788f3fbd535SBarry Smith   }
789f3fbd535SBarry Smith 
790f3fbd535SBarry Smith   if (viewer) {
791f3fbd535SBarry Smith     for (i=1; i<n; i++) {
792f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
793f3fbd535SBarry Smith     }
794f3fbd535SBarry Smith     for (i=0; i<n; i++) {
795f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
796f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
797f3fbd535SBarry Smith     }
798f3fbd535SBarry Smith   }
799f3fbd535SBarry Smith   PetscFunctionReturn(0);
800f3fbd535SBarry Smith }
801f3fbd535SBarry Smith 
802f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
803f3fbd535SBarry Smith 
804f3fbd535SBarry Smith #undef __FUNCT__
8059dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels"
8064b9ad928SBarry Smith /*@
80797177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
8084b9ad928SBarry Smith 
8094b9ad928SBarry Smith    Not Collective
8104b9ad928SBarry Smith 
8114b9ad928SBarry Smith    Input Parameter:
8124b9ad928SBarry Smith .  pc - the preconditioner context
8134b9ad928SBarry Smith 
8144b9ad928SBarry Smith    Output parameter:
8154b9ad928SBarry Smith .  levels - the number of levels
8164b9ad928SBarry Smith 
8174b9ad928SBarry Smith    Level: advanced
8184b9ad928SBarry Smith 
8194b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8204b9ad928SBarry Smith 
82197177400SBarry Smith .seealso: PCMGSetLevels()
8224b9ad928SBarry Smith @*/
8237087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8244b9ad928SBarry Smith {
825f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8264b9ad928SBarry Smith 
8274b9ad928SBarry Smith   PetscFunctionBegin;
8280700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8294482741eSBarry Smith   PetscValidIntPointer(levels,2);
830f3fbd535SBarry Smith   *levels = mg->nlevels;
8314b9ad928SBarry Smith   PetscFunctionReturn(0);
8324b9ad928SBarry Smith }
8334b9ad928SBarry Smith 
8344b9ad928SBarry Smith #undef __FUNCT__
8359dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType"
8364b9ad928SBarry Smith /*@
83797177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8384b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8394b9ad928SBarry Smith 
840ad4df100SBarry Smith    Logically Collective on PC
8414b9ad928SBarry Smith 
8424b9ad928SBarry Smith    Input Parameters:
8434b9ad928SBarry Smith +  pc - the preconditioner context
8449dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
8459dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
8464b9ad928SBarry Smith 
8474b9ad928SBarry Smith    Options Database Key:
8484b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
8494b9ad928SBarry Smith    additive, full, kaskade
8504b9ad928SBarry Smith 
8514b9ad928SBarry Smith    Level: advanced
8524b9ad928SBarry Smith 
8534b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
8544b9ad928SBarry Smith 
85597177400SBarry Smith .seealso: PCMGSetLevels()
8564b9ad928SBarry Smith @*/
8577087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
8584b9ad928SBarry Smith {
859f3fbd535SBarry Smith   PC_MG                   *mg = (PC_MG*)pc->data;
8604b9ad928SBarry Smith 
8614b9ad928SBarry Smith   PetscFunctionBegin;
8620700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
863c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
86431567311SBarry Smith   mg->am = form;
8659dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
8664b9ad928SBarry Smith   else pc->ops->applyrichardson = 0;
8674b9ad928SBarry Smith   PetscFunctionReturn(0);
8684b9ad928SBarry Smith }
8694b9ad928SBarry Smith 
8704b9ad928SBarry Smith #undef __FUNCT__
8710d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType"
8724b9ad928SBarry Smith /*@
8730d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
8744b9ad928SBarry Smith    complicated cycling.
8754b9ad928SBarry Smith 
876ad4df100SBarry Smith    Logically Collective on PC
8774b9ad928SBarry Smith 
8784b9ad928SBarry Smith    Input Parameters:
879c2be2410SBarry Smith +  pc - the multigrid context
8800d353602SBarry Smith -  PC_MG_CYCLE_V or PC_MG_CYCLE_W
8814b9ad928SBarry Smith 
8824b9ad928SBarry Smith    Options Database Key:
8830d353602SBarry Smith $  -pc_mg_cycle_type v or w
8844b9ad928SBarry Smith 
8854b9ad928SBarry Smith    Level: advanced
8864b9ad928SBarry Smith 
8874b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
8884b9ad928SBarry Smith 
8890d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
8904b9ad928SBarry Smith @*/
8917087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
8924b9ad928SBarry Smith {
893f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
894f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
89579416396SBarry Smith   PetscInt     i,levels;
8964b9ad928SBarry Smith 
8974b9ad928SBarry Smith   PetscFunctionBegin;
8980700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
899e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
900c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
901f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9024b9ad928SBarry Smith 
9034b9ad928SBarry Smith   for (i=0; i<levels; i++) {
904f3fbd535SBarry Smith     mglevels[i]->cycles  = n;
9054b9ad928SBarry Smith   }
9064b9ad928SBarry Smith   PetscFunctionReturn(0);
9074b9ad928SBarry Smith }
9084b9ad928SBarry Smith 
9094b9ad928SBarry Smith #undef __FUNCT__
9108cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles"
9118cc2d5dfSBarry Smith /*@
9128cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
9138cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9148cc2d5dfSBarry Smith 
915ad4df100SBarry Smith    Logically Collective on PC
9168cc2d5dfSBarry Smith 
9178cc2d5dfSBarry Smith    Input Parameters:
9188cc2d5dfSBarry Smith +  pc - the multigrid context
9198cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
9208cc2d5dfSBarry Smith 
9218cc2d5dfSBarry Smith    Options Database Key:
9228cc2d5dfSBarry Smith $  -pc_mg_multiplicative_cycles n
9238cc2d5dfSBarry Smith 
9248cc2d5dfSBarry Smith    Level: advanced
9258cc2d5dfSBarry Smith 
9268cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
9278cc2d5dfSBarry Smith 
9288cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9298cc2d5dfSBarry Smith 
9308cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
9318cc2d5dfSBarry Smith @*/
9327087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
9338cc2d5dfSBarry Smith {
934f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
935f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
9368cc2d5dfSBarry Smith   PetscInt     i,levels;
9378cc2d5dfSBarry Smith 
9388cc2d5dfSBarry Smith   PetscFunctionBegin;
9390700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
940e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
941c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
942f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9438cc2d5dfSBarry Smith 
9448cc2d5dfSBarry Smith   for (i=0; i<levels; i++) {
94531567311SBarry Smith     mg->cyclesperpcapply  = n;
9468cc2d5dfSBarry Smith   }
9478cc2d5dfSBarry Smith   PetscFunctionReturn(0);
9488cc2d5dfSBarry Smith }
9498cc2d5dfSBarry Smith 
9508cc2d5dfSBarry Smith #undef __FUNCT__
9519dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin"
952c2be2410SBarry Smith /*@
95397177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
954c2be2410SBarry Smith       finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
955c2be2410SBarry Smith 
956ad4df100SBarry Smith    Logically Collective on PC
957c2be2410SBarry Smith 
958c2be2410SBarry Smith    Input Parameters:
959c91913d3SJed Brown +  pc - the multigrid context
960c91913d3SJed Brown -  use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
961c2be2410SBarry Smith 
962c2be2410SBarry Smith    Options Database Key:
963c2be2410SBarry Smith $  -pc_mg_galerkin
964c2be2410SBarry Smith 
965c2be2410SBarry Smith    Level: intermediate
966c2be2410SBarry Smith 
967c2be2410SBarry Smith .keywords: MG, set, Galerkin
968c2be2410SBarry Smith 
9693fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin()
9703fc8bf9cSBarry Smith 
971c2be2410SBarry Smith @*/
972c91913d3SJed Brown PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
973c2be2410SBarry Smith {
974f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
975c2be2410SBarry Smith 
976c2be2410SBarry Smith   PetscFunctionBegin;
9770700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
978789726d7SBarry Smith   mg->galerkin = use ? 1 : 0;
979c2be2410SBarry Smith   PetscFunctionReturn(0);
980c2be2410SBarry Smith }
981c2be2410SBarry Smith 
982c2be2410SBarry Smith #undef __FUNCT__
9833fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin"
9843fc8bf9cSBarry Smith /*@
9853fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
9863fc8bf9cSBarry Smith       A_i-1 = r_i * A_i * r_i^t
9873fc8bf9cSBarry Smith 
9883fc8bf9cSBarry Smith    Not Collective
9893fc8bf9cSBarry Smith 
9903fc8bf9cSBarry Smith    Input Parameter:
9913fc8bf9cSBarry Smith .  pc - the multigrid context
9923fc8bf9cSBarry Smith 
9933fc8bf9cSBarry Smith    Output Parameter:
9943fc8bf9cSBarry Smith .  gelerkin - PETSC_TRUE or PETSC_FALSE
9953fc8bf9cSBarry Smith 
9963fc8bf9cSBarry Smith    Options Database Key:
9973fc8bf9cSBarry Smith $  -pc_mg_galerkin
9983fc8bf9cSBarry Smith 
9993fc8bf9cSBarry Smith    Level: intermediate
10003fc8bf9cSBarry Smith 
10013fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
10023fc8bf9cSBarry Smith 
10033fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin()
10043fc8bf9cSBarry Smith 
10053fc8bf9cSBarry Smith @*/
10067087cfbeSBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PetscBool  *galerkin)
10073fc8bf9cSBarry Smith {
1008f3fbd535SBarry Smith   PC_MG        *mg = (PC_MG*)pc->data;
10093fc8bf9cSBarry Smith 
10103fc8bf9cSBarry Smith   PetscFunctionBegin;
10110700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
101213fdb3acSJose Roman   *galerkin = (PetscBool)mg->galerkin;
10133fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10143fc8bf9cSBarry Smith }
10153fc8bf9cSBarry Smith 
10163fc8bf9cSBarry Smith #undef __FUNCT__
10179dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown"
10184b9ad928SBarry Smith /*@
101997177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
102097177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
10214b9ad928SBarry Smith    pre-smoothing steps on different levels.
10224b9ad928SBarry Smith 
1023ad4df100SBarry Smith    Logically Collective on PC
10244b9ad928SBarry Smith 
10254b9ad928SBarry Smith    Input Parameters:
10264b9ad928SBarry Smith +  mg - the multigrid context
10274b9ad928SBarry Smith -  n - the number of smoothing steps
10284b9ad928SBarry Smith 
10294b9ad928SBarry Smith    Options Database Key:
10304b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
10314b9ad928SBarry Smith 
10324b9ad928SBarry Smith    Level: advanced
10334b9ad928SBarry Smith 
10344b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
10354b9ad928SBarry Smith 
103697177400SBarry Smith .seealso: PCMGSetNumberSmoothUp()
10374b9ad928SBarry Smith @*/
10387087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
10394b9ad928SBarry Smith {
1040f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1041f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10426849ba73SBarry Smith   PetscErrorCode ierr;
104379416396SBarry Smith   PetscInt       i,levels;
10444b9ad928SBarry Smith 
10454b9ad928SBarry Smith   PetscFunctionBegin;
10460700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1047e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1048c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1049f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10504b9ad928SBarry Smith 
1051b05257ddSBarry Smith   for (i=1; i<levels; i++) {
10524b9ad928SBarry Smith     /* make sure smoother up and down are different */
105397177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1054f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
105531567311SBarry Smith     mg->default_smoothd = n;
10564b9ad928SBarry Smith   }
10574b9ad928SBarry Smith   PetscFunctionReturn(0);
10584b9ad928SBarry Smith }
10594b9ad928SBarry Smith 
10604b9ad928SBarry Smith #undef __FUNCT__
10619dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp"
10624b9ad928SBarry Smith /*@
106397177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
106497177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
10654b9ad928SBarry Smith    post-smoothing steps on different levels.
10664b9ad928SBarry Smith 
1067ad4df100SBarry Smith    Logically Collective on PC
10684b9ad928SBarry Smith 
10694b9ad928SBarry Smith    Input Parameters:
10704b9ad928SBarry Smith +  mg - the multigrid context
10714b9ad928SBarry Smith -  n - the number of smoothing steps
10724b9ad928SBarry Smith 
10734b9ad928SBarry Smith    Options Database Key:
10744b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
10754b9ad928SBarry Smith 
10764b9ad928SBarry Smith    Level: advanced
10774b9ad928SBarry Smith 
10784b9ad928SBarry Smith    Note: this does not set a value on the coarsest grid, since we assume that
1079a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
10804b9ad928SBarry Smith 
10814b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
10824b9ad928SBarry Smith 
108397177400SBarry Smith .seealso: PCMGSetNumberSmoothDown()
10844b9ad928SBarry Smith @*/
10857087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
10864b9ad928SBarry Smith {
1087f3fbd535SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1088f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
10896849ba73SBarry Smith   PetscErrorCode ierr;
109079416396SBarry Smith   PetscInt       i,levels;
10914b9ad928SBarry Smith 
10924b9ad928SBarry Smith   PetscFunctionBegin;
10930700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1094e7e72b3dSBarry Smith   if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1095c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1096f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10974b9ad928SBarry Smith 
10984b9ad928SBarry Smith   for (i=1; i<levels; i++) {
10994b9ad928SBarry Smith     /* make sure smoother up and down are different */
110097177400SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr);
1101f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
110231567311SBarry Smith     mg->default_smoothu = n;
11034b9ad928SBarry Smith   }
11044b9ad928SBarry Smith   PetscFunctionReturn(0);
11054b9ad928SBarry Smith }
11064b9ad928SBarry Smith 
11074b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
11084b9ad928SBarry Smith 
11093b09bd56SBarry Smith /*MC
1110ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
11113b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
11123b09bd56SBarry Smith 
11133b09bd56SBarry Smith    Options Database Keys:
11143b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
11150d353602SBarry Smith .  -pc_mg_cycles v or w
111679416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
11173b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
11188c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
11193b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
11203b09bd56SBarry Smith .  -pc_mg_monitor - print information on the multigrid convergence
11218cf6037eSBarry Smith .  -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
11228cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
11238cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1124e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
11258cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
11268cf6037eSBarry Smith                         to the binary output file called binaryoutput
11273b09bd56SBarry Smith 
112824c3aa18SBarry Smith    Notes: By default this uses GMRES on the fine grid smoother so this should be used with KSPFGMRES or the smoother changed to not use GMRES
11293b09bd56SBarry Smith 
11308cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
11318cf6037eSBarry Smith 
11323b09bd56SBarry Smith    Level: intermediate
11333b09bd56SBarry Smith 
11348f87f92bSBarry Smith    Concepts: multigrid/multilevel
11353b09bd56SBarry Smith 
11368cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
11370d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
113897177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
113997177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
11400d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
11413b09bd56SBarry Smith M*/
11423b09bd56SBarry Smith 
11434b9ad928SBarry Smith EXTERN_C_BEGIN
11444b9ad928SBarry Smith #undef __FUNCT__
11454b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG"
11467087cfbeSBarry Smith PetscErrorCode  PCCreate_MG(PC pc)
11474b9ad928SBarry Smith {
1148f3fbd535SBarry Smith   PC_MG          *mg;
1149f3fbd535SBarry Smith   PetscErrorCode ierr;
1150f3fbd535SBarry Smith 
11514b9ad928SBarry Smith   PetscFunctionBegin;
1152f3fbd535SBarry Smith   ierr        = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr);
1153f3fbd535SBarry Smith   pc->data    = (void*)mg;
1154f3fbd535SBarry Smith   mg->nlevels = -1;
1155f3fbd535SBarry Smith 
11564b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
11574b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1158a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
11594b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
11604b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
11614b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
11624b9ad928SBarry Smith   PetscFunctionReturn(0);
11634b9ad928SBarry Smith }
11644b9ad928SBarry Smith EXTERN_C_END
1165