xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision c0decd05c6848b80907752eef350b55c8c90e696)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>                    /*I "petscksp.h" I*/
61e25c274SJed Brown #include <petscdm.h>
7391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*);
84b9ad928SBarry Smith 
931567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
104b9ad928SBarry Smith {
1131567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1231567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
136849ba73SBarry Smith   PetscErrorCode ierr;
1431567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
154b9ad928SBarry Smith 
164b9ad928SBarry Smith   PetscFunctionBegin;
1763e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
1831567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
19*c0decd05SBarry Smith   ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr);
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;
3357622a8eSBarry Smith           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr);
344b9ad928SBarry Smith         } else {
354d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
3657622a8eSBarry Smith           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)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 */
55*c0decd05SBarry Smith     ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr);
5663e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith   }
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
61ace3abfcSBarry 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)
624b9ad928SBarry Smith {
63f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
64f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
65dfbe8321SBarry Smith   PetscErrorCode ierr;
66391689abSStefano Zampini   PC             tpc;
67391689abSStefano Zampini   PetscBool      changeu,changed;
68f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
694b9ad928SBarry Smith 
704b9ad928SBarry Smith   PetscFunctionBegin;
71a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
72a762d673SBarry Smith   for (i=0; i<levels; i++) {
73a762d673SBarry Smith     if (!mglevels[i]->A) {
74a762d673SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
75a762d673SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
76a762d673SBarry Smith     }
77a762d673SBarry Smith   }
78391689abSStefano Zampini 
79391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
80391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
81391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
82391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
83391689abSStefano Zampini   if (!changed && !changeu) {
84391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
85f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
86391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
87391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
88391689abSStefano Zampini       Vec *vec;
89391689abSStefano Zampini 
90391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
91391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
927ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
93391689abSStefano Zampini     }
94391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
95391689abSStefano Zampini   }
96f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
974b9ad928SBarry Smith 
9831567311SBarry Smith   mg->rtol   = rtol;
9931567311SBarry Smith   mg->abstol = abstol;
10031567311SBarry Smith   mg->dtol   = dtol;
1014b9ad928SBarry Smith   if (rtol) {
1024b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
1034b9ad928SBarry Smith     PetscReal rnorm;
1047319c654SBarry Smith     if (zeroguess) {
1057319c654SBarry Smith       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
1067319c654SBarry Smith     } else {
107f3fbd535SBarry Smith       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
1084b9ad928SBarry Smith       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
1097319c654SBarry Smith     }
11031567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
1112fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
1122fa5cd67SKarl Rupp   else mg->ttol = 0.0;
1134b9ad928SBarry Smith 
1144d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
115336babb1SJed Brown      stop prematurely due to small residual */
1164d0a8057SBarry Smith   for (i=1; i<levels; i++) {
117f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
118f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
11923067569SBarry Smith       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
12023067569SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
121f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1224b9ad928SBarry Smith     }
1234d0a8057SBarry Smith   }
1244d0a8057SBarry Smith 
1254d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1264d0a8057SBarry Smith   for (i=0; i<its; i++) {
127f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1284d0a8057SBarry Smith     if (*reason) break;
1294d0a8057SBarry Smith   }
1304d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1314d0a8057SBarry Smith   *outits = i;
132391689abSStefano Zampini   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
1334b9ad928SBarry Smith   PetscFunctionReturn(0);
1344b9ad928SBarry Smith }
1354b9ad928SBarry Smith 
1363aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1373aeaf226SBarry Smith {
1383aeaf226SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1393aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1403aeaf226SBarry Smith   PetscErrorCode ierr;
1413aeaf226SBarry Smith   PetscInt       i,n;
1423aeaf226SBarry Smith 
1433aeaf226SBarry Smith   PetscFunctionBegin;
1443aeaf226SBarry Smith   if (mglevels) {
1453aeaf226SBarry Smith     n = mglevels[0]->levels;
1463aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1473aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1483aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1493aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
1503aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1513aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
1520de8493bSLawrence Mitchell       ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr);
15373250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1543aeaf226SBarry Smith     }
155391689abSStefano Zampini     /* this is not null only if the smoother on the finest level
156391689abSStefano Zampini        changes the rhs during PreSolve */
157391689abSStefano Zampini     ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr);
1583aeaf226SBarry Smith 
1593aeaf226SBarry Smith     for (i=0; i<n; i++) {
1603aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1613aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1623aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1633aeaf226SBarry Smith       }
1643aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1653aeaf226SBarry Smith     }
1663aeaf226SBarry Smith   }
1673aeaf226SBarry Smith   PetscFunctionReturn(0);
1683aeaf226SBarry Smith }
1693aeaf226SBarry Smith 
1704b9ad928SBarry Smith /*@C
17197177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1724b9ad928SBarry Smith    Must be called before any other MG routine.
1734b9ad928SBarry Smith 
174ad4df100SBarry Smith    Logically Collective on PC
1754b9ad928SBarry Smith 
1764b9ad928SBarry Smith    Input Parameters:
1774b9ad928SBarry Smith +  pc - the preconditioner context
1784b9ad928SBarry Smith .  levels - the number of levels
1794b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1802bf68e3eSBarry Smith            on smaller sets of processors.
1814b9ad928SBarry Smith 
1824b9ad928SBarry Smith    Level: intermediate
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith    Notes:
1854b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1864b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1874b9ad928SBarry Smith 
1884b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1894b9ad928SBarry Smith 
19097177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1914b9ad928SBarry Smith @*/
1927087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1934b9ad928SBarry Smith {
194dfbe8321SBarry Smith   PetscErrorCode ierr;
195f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
196ce94432eSBarry Smith   MPI_Comm       comm;
1973aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
19810eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
19910167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
200f3fbd535SBarry Smith   PetscInt       i;
201f3fbd535SBarry Smith   PetscMPIInt    size;
202f3fbd535SBarry Smith   const char     *prefix;
203f3fbd535SBarry Smith   PC             ipc;
2043aeaf226SBarry Smith   PetscInt       n;
2054b9ad928SBarry Smith 
2064b9ad928SBarry Smith   PetscFunctionBegin;
2070700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
208548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
209548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
2101a97d7b6SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
2113aeaf226SBarry Smith   if (mglevels) {
21210eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
2133aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
2143aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
2153aeaf226SBarry Smith     n    = mglevels[0]->levels;
2163aeaf226SBarry Smith     for (i=0; i<n; i++) {
2173aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2183aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
2193aeaf226SBarry Smith       }
2203aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
2213aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
2223aeaf226SBarry Smith     }
2233aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
2243aeaf226SBarry Smith   }
225f3fbd535SBarry Smith 
226f3fbd535SBarry Smith   mg->nlevels = levels;
227f3fbd535SBarry Smith 
228785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
2293bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
230f3fbd535SBarry Smith 
231f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
232f3fbd535SBarry Smith 
233a9db87a2SMatthew G Knepley   mg->stageApply = 0;
234f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
235b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
2362fa5cd67SKarl Rupp 
237f3fbd535SBarry Smith     mglevels[i]->level               = i;
238f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
23910eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
240336babb1SJed Brown     mg->default_smoothu              = 2;
241336babb1SJed Brown     mg->default_smoothd              = 2;
24263e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
24363e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
24463e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
24563e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
246f3fbd535SBarry Smith 
247f3fbd535SBarry Smith     if (comms) comm = comms[i];
248f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
249422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
2500ee9a668SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
2510ee9a668SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
2520ee9a668SBarry Smith     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
2530ee9a668SBarry Smith     if (i || levels == 1) {
2540ee9a668SBarry Smith       char tprefix[128];
2550ee9a668SBarry Smith 
256336babb1SJed Brown       ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2570059c7bdSJed Brown       ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
258336babb1SJed Brown       ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
259336babb1SJed Brown       ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
260336babb1SJed Brown       ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
2610ee9a668SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
262f3fbd535SBarry Smith 
2630ee9a668SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
2640ee9a668SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
2650ee9a668SBarry Smith     } else {
266f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
267f3fbd535SBarry Smith 
2689dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
269f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
270f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
271f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
272f3fbd535SBarry Smith       if (size > 1) {
273f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
274f3fbd535SBarry Smith       } else {
275f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
276f3fbd535SBarry Smith       }
277753b7fb9SBarry Smith       ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
278f3fbd535SBarry Smith     }
2793bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
2802fa5cd67SKarl Rupp 
281f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
28231567311SBarry Smith     mg->rtol             = 0.0;
28331567311SBarry Smith     mg->abstol           = 0.0;
28431567311SBarry Smith     mg->dtol             = 0.0;
28531567311SBarry Smith     mg->ttol             = 0.0;
28631567311SBarry Smith     mg->cyclesperpcapply = 1;
287f3fbd535SBarry Smith   }
288f3fbd535SBarry Smith   mg->levels = mglevels;
28910eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2904b9ad928SBarry Smith   PetscFunctionReturn(0);
2914b9ad928SBarry Smith }
2924b9ad928SBarry Smith 
293c07bf074SBarry Smith 
294c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
295c07bf074SBarry Smith {
296c07bf074SBarry Smith   PetscErrorCode ierr;
297a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
298a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
299a06653b4SBarry Smith   PetscInt       i,n;
300c07bf074SBarry Smith 
301c07bf074SBarry Smith   PetscFunctionBegin;
302a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
303a06653b4SBarry Smith   if (mglevels) {
304a06653b4SBarry Smith     n = mglevels[0]->levels;
305a06653b4SBarry Smith     for (i=0; i<n; i++) {
306a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
3076bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
308a06653b4SBarry Smith       }
3096bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
310a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
311a06653b4SBarry Smith     }
312a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
313a06653b4SBarry Smith   }
314c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
315f3fbd535SBarry Smith   PetscFunctionReturn(0);
316f3fbd535SBarry Smith }
317f3fbd535SBarry Smith 
318f3fbd535SBarry Smith 
319f3fbd535SBarry Smith 
32009573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
32109573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
32209573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
323f3fbd535SBarry Smith 
324f3fbd535SBarry Smith /*
325f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
326f3fbd535SBarry Smith              or full cycle of multigrid.
327f3fbd535SBarry Smith 
328f3fbd535SBarry Smith   Note:
329f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
330f3fbd535SBarry Smith */
331f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
332f3fbd535SBarry Smith {
333f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
334f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
335f3fbd535SBarry Smith   PetscErrorCode ierr;
336391689abSStefano Zampini   PC             tpc;
337f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
338391689abSStefano Zampini   PetscBool      changeu,changed;
339f3fbd535SBarry Smith 
340f3fbd535SBarry Smith   PetscFunctionBegin;
341a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
342e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
343298cc208SBarry Smith   for (i=0; i<levels; i++) {
344e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
34523ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
346298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
347e1d8e5deSBarry Smith     }
348e1d8e5deSBarry Smith   }
349e1d8e5deSBarry Smith 
350391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
351391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
352391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
353391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
354391689abSStefano Zampini   if (!changeu && !changed) {
355391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
356f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
357391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
358391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
359391689abSStefano Zampini       Vec *vec;
360391689abSStefano Zampini 
361391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
362391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
3637ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
364391689abSStefano Zampini     }
365391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
366391689abSStefano Zampini   }
367f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
368391689abSStefano Zampini 
36931567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
370f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
37131567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
3720298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
373f3fbd535SBarry Smith     }
3742fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
37531567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
3762fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
37731567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
3782fa5cd67SKarl Rupp   } else {
379f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
380f3fbd535SBarry Smith   }
381a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
382391689abSStefano Zampini   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
383f3fbd535SBarry Smith   PetscFunctionReturn(0);
384f3fbd535SBarry Smith }
385f3fbd535SBarry Smith 
386f3fbd535SBarry Smith 
3874416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
388f3fbd535SBarry Smith {
389f3fbd535SBarry Smith   PetscErrorCode   ierr;
390710315b6SLawrence Mitchell   PetscInt         levels,cycles;
3912134b1e4SBarry Smith   PetscBool        flg;
392f3fbd535SBarry Smith   PC_MG            *mg = (PC_MG*)pc->data;
3933d3eaba7SBarry Smith   PC_MG_Levels     **mglevels;
394f3fbd535SBarry Smith   PCMGType         mgtype;
395f3fbd535SBarry Smith   PCMGCycleType    mgctype;
3962134b1e4SBarry Smith   PCMGGalerkinType gtype;
397f3fbd535SBarry Smith 
398f3fbd535SBarry Smith   PetscFunctionBegin;
3991d743356SStefano Zampini   levels = PetscMax(mg->nlevels,1);
400e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
401f3fbd535SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
4021a97d7b6SStefano Zampini   if (!flg && !mg->levels && pc->dm) {
403eb3f98d2SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
404eb3f98d2SBarry Smith     levels++;
4053aeaf226SBarry Smith     mg->usedmfornumberoflevels = PETSC_TRUE;
406eb3f98d2SBarry Smith   }
4070298fd71SBarry Smith   ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
4083aeaf226SBarry Smith   mglevels = mg->levels;
4093aeaf226SBarry Smith 
410f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
411f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
412f3fbd535SBarry Smith   if (flg) {
413f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
4142fa5cd67SKarl Rupp   }
4152134b1e4SBarry Smith   gtype = mg->galerkin;
4162134b1e4SBarry Smith   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
4172134b1e4SBarry Smith   if (flg) {
4182134b1e4SBarry Smith     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
419f3fbd535SBarry Smith   }
420f56b1265SBarry Smith   flg = PETSC_FALSE;
421f442ab6aSBarry Smith   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
422f442ab6aSBarry Smith   if (flg) {
423f442ab6aSBarry Smith     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
424f442ab6aSBarry Smith   }
42531567311SBarry Smith   mgtype = mg->am;
426f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
427f3fbd535SBarry Smith   if (flg) {
428f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
429f3fbd535SBarry Smith   }
43031567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
4315363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
432f3fbd535SBarry Smith     if (flg) {
433f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
434f3fbd535SBarry Smith     }
435f3fbd535SBarry Smith   }
436f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4370298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
438f3fbd535SBarry Smith   if (flg) {
439f3fbd535SBarry Smith     PetscInt i;
440f3fbd535SBarry Smith     char     eventname[128];
4411a97d7b6SStefano Zampini 
442f3fbd535SBarry Smith     levels = mglevels[0]->levels;
443f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
444f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
44563e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
446f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
44763e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
448f3fbd535SBarry Smith       if (i) {
449f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
45063e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
451f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
45263e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
453f3fbd535SBarry Smith       }
454f3fbd535SBarry Smith     }
455eec35531SMatthew G Knepley 
456ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
457eec35531SMatthew G Knepley     {
458eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
459eec35531SMatthew G Knepley       PetscStageLog stageLog;
460eec35531SMatthew G Knepley       PetscInt      st;
461eec35531SMatthew G Knepley 
462eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
463eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
464eec35531SMatthew G Knepley         PetscBool same;
465eec35531SMatthew G Knepley 
466eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4672fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
468eec35531SMatthew G Knepley       }
469eec35531SMatthew G Knepley       if (!mg->stageApply) {
470eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
471eec35531SMatthew G Knepley       }
472eec35531SMatthew G Knepley     }
473ec60ca73SSatish Balay #endif
474f3fbd535SBarry Smith   }
475f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
476f3fbd535SBarry Smith   PetscFunctionReturn(0);
477f3fbd535SBarry Smith }
478f3fbd535SBarry Smith 
4796a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4806a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
4812134b1e4SBarry Smith const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0};
482f3fbd535SBarry Smith 
4839804daf3SBarry Smith #include <petscdraw.h>
484f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
485f3fbd535SBarry Smith {
486f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
487f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
488f3fbd535SBarry Smith   PetscErrorCode ierr;
489e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
490219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
491f3fbd535SBarry Smith 
492f3fbd535SBarry Smith   PetscFunctionBegin;
493251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4945b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
495219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
496f3fbd535SBarry Smith   if (iascii) {
497e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
498efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
49931567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
50031567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
501f3fbd535SBarry Smith     }
5022134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
503f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
5042134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
5052134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
5062134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
5072134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
5082134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
5092134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
5104f66f45eSBarry Smith     } else {
5114f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
512f3fbd535SBarry Smith     }
5135adeb434SBarry Smith     if (mg->view){
5145adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
5155adeb434SBarry Smith     }
516f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
517f3fbd535SBarry Smith       if (!i) {
51889cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
519f3fbd535SBarry Smith       } else {
52089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
521f3fbd535SBarry Smith       }
522f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
523f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
524f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
525f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
526f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
527f3fbd535SBarry Smith       } else if (i) {
52889cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
529f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
530f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
531f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
532f3fbd535SBarry Smith       }
533f3fbd535SBarry Smith     }
5345b0b0462SBarry Smith   } else if (isbinary) {
5355b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5365b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5375b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5385b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5395b0b0462SBarry Smith       }
5405b0b0462SBarry Smith     }
541219b1045SBarry Smith   } else if (isdraw) {
542219b1045SBarry Smith     PetscDraw draw;
543b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
544219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
545219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5460298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
547219b1045SBarry Smith     bottom = y - th;
548219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
549b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
550219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
551219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
552219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
553b4375e8dSPeter Brune       } else {
554b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
555b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
556b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
557b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
558b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
559b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
560b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
561b4375e8dSPeter Brune       }
5620298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5631cd381d6SBarry Smith       bottom -= th;
564219b1045SBarry Smith     }
5655b0b0462SBarry Smith   }
566f3fbd535SBarry Smith   PetscFunctionReturn(0);
567f3fbd535SBarry Smith }
568f3fbd535SBarry Smith 
569af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
570af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
571cab2e9ccSBarry Smith 
572f3fbd535SBarry Smith /*
573f3fbd535SBarry Smith     Calls setup for the KSP on each level
574f3fbd535SBarry Smith */
575f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
576f3fbd535SBarry Smith {
577f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
578f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
579f3fbd535SBarry Smith   PetscErrorCode ierr;
5801a97d7b6SStefano Zampini   PetscInt       i,n;
58198e04984SBarry Smith   PC             cpc;
5823db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
583f3fbd535SBarry Smith   Mat            dA,dB;
584f3fbd535SBarry Smith   Vec            tvec;
585218a07d4SBarry Smith   DM             *dms;
586649052a6SBarry Smith   PetscViewer    viewer = 0;
5872134b1e4SBarry Smith   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
588f3fbd535SBarry Smith 
589f3fbd535SBarry Smith   PetscFunctionBegin;
5901a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
5911a97d7b6SStefano Zampini   n = mglevels[0]->levels;
59201bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5933aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5943aeaf226SBarry Smith     PetscInt levels;
5953aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5963aeaf226SBarry Smith     levels++;
5973aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5980298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5993aeaf226SBarry Smith       n        = levels;
6003aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
6013aeaf226SBarry Smith       mglevels =  mg->levels;
6023aeaf226SBarry Smith     }
6033aeaf226SBarry Smith   }
60498e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
6053aeaf226SBarry Smith 
606f3fbd535SBarry Smith 
607f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
608f3fbd535SBarry Smith   /* so use those from global PC */
609f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
6100298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
61198e04984SBarry Smith   if (opsset) {
61298e04984SBarry Smith     Mat mmat;
61323ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
61498e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
61598e04984SBarry Smith   }
6164a5f07a5SMark F. Adams 
6174a5f07a5SMark F. Adams   if (!opsset) {
61871b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
6194a5f07a5SMark F. Adams     if(use_amat){
620f3fbd535SBarry 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);
62123ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
622f3fbd535SBarry Smith     }
6234a5f07a5SMark F. Adams     else {
6244a5f07a5SMark F. Adams       ierr = PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr);
62523ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
6264a5f07a5SMark F. Adams     }
6274a5f07a5SMark F. Adams   }
628f3fbd535SBarry Smith 
6299c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
6309c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
6319c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
6329c7ffb39SBarry Smith       continue;
6339c7ffb39SBarry Smith     }
6349c7ffb39SBarry Smith   }
6352134b1e4SBarry Smith 
6362134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
6372134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
6382134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
6392134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
6402134b1e4SBarry Smith   }
6412134b1e4SBarry Smith 
6422134b1e4SBarry Smith 
6439c7ffb39SBarry Smith   /*
6449c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6459c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6469c7ffb39SBarry Smith   */
6472134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
6482d2b81a6SBarry Smith 	/* construct the interpolation from the DMs */
6492e499ae9SBarry Smith     Mat p;
65073250ac0SBarry Smith     Vec rscale;
651785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
652218a07d4SBarry Smith     dms[n-1] = pc->dm;
653ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
654ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
65541fce8fdSHong Zhang 	/*
65641fce8fdSHong Zhang 	   Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
65741fce8fdSHong Zhang 	   Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
65841fce8fdSHong Zhang 	   But it is safe to use -dm_mat_type.
65941fce8fdSHong Zhang 
66041fce8fdSHong Zhang 	   The mat type should not be hardcoded like this, we need to find a better way.
66141fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
66241fce8fdSHong Zhang     */
663218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
664942e3340SBarry Smith       DMKSP     kdm;
665eab52d2bSLawrence Mitchell       PetscBool dmhasrestrict, dmhasinject;
6663c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6672134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
668942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
669d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
670d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6710298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6720298fd71SBarry Smith       kdm->rhsctx          = NULL;
67324c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
674e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6752d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
67600ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
67773250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6786bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
679218a07d4SBarry Smith       }
6803ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
6813ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
6823ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
6833ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
6843ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
6853ad4599aSBarry Smith       }
686eab52d2bSLawrence Mitchell       ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr);
687eab52d2bSLawrence Mitchell       if (dmhasinject && !mglevels[i+1]->inject){
688eab52d2bSLawrence Mitchell         ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr);
689eab52d2bSLawrence Mitchell         ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr);
690eab52d2bSLawrence Mitchell         ierr = MatDestroy(&p);CHKERRQ(ierr);
691eab52d2bSLawrence Mitchell       }
69224c3aa18SBarry Smith     }
6932d2b81a6SBarry Smith 
694ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
6952d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
696ce4cda84SJed Brown   }
697cab2e9ccSBarry Smith 
698ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
6992134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
700cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
701cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
702218a07d4SBarry Smith   }
703218a07d4SBarry Smith 
7042134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
7052134b1e4SBarry Smith     Mat       A,B;
7062134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
7072134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
7082134b1e4SBarry Smith 
7092134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
7102134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
7112134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
712f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
713b9367914SBarry Smith       if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0");
714b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
715b9367914SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
716b9367914SBarry Smith       }
717b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
718b9367914SBarry Smith         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
719b9367914SBarry Smith       }
7202134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
7212134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
7222134b1e4SBarry Smith       }
7232134b1e4SBarry Smith       if (doA) {
7242df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
7252134b1e4SBarry Smith       }
7262134b1e4SBarry Smith       if (doB) {
7272df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
728b9367914SBarry Smith       }
7292134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
7302134b1e4SBarry Smith       if (!doA && dAeqdB) {
7312134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
7322134b1e4SBarry Smith         A = B;
7332134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX ) {
7342134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
7352134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
736b9367914SBarry Smith       }
7372134b1e4SBarry Smith       if (!doB && dAeqdB) {
7382134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
7392134b1e4SBarry Smith         B = A;
7402134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
74123ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
7422134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
7437d537d92SJed Brown       }
7442134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
7452134b1e4SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
7462134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
7472134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
7482134b1e4SBarry Smith       }
7492134b1e4SBarry Smith       dA = A;
750f3fbd535SBarry Smith       dB = B;
751f3fbd535SBarry Smith     }
752f3fbd535SBarry Smith   }
7532134b1e4SBarry Smith   if (needRestricts && pc->dm && pc->dm->x) {
754cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
755cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
756c88c5224SJed Brown       Mat R;
757c88c5224SJed Brown       Vec rscale;
758cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
759cab2e9ccSBarry Smith         Vec *vecs;
7602a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
761cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
762cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
763cab2e9ccSBarry Smith       }
764c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
765c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
766c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
767c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
768cab2e9ccSBarry Smith     }
769f3fbd535SBarry Smith   }
7702134b1e4SBarry Smith   if (needRestricts && pc->dm) {
771caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
772ccceb50cSJed Brown       DM  dmfine,dmcoarse;
773caa4e7f2SJed Brown       Mat Restrict,Inject;
774caa4e7f2SJed Brown       Vec rscale;
775ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
776ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
777caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
778caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
779eab52d2bSLawrence Mitchell       ierr   = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr);
780caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
781caa4e7f2SJed Brown     }
782caa4e7f2SJed Brown   }
783f3fbd535SBarry Smith 
784f3fbd535SBarry Smith   if (!pc->setupcalled) {
785f3fbd535SBarry Smith     for (i=0; i<n; i++) {
786f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
787f3fbd535SBarry Smith     }
788f3fbd535SBarry Smith     for (i=1; i<n; i++) {
789f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
790f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
791f3fbd535SBarry Smith       }
792f3fbd535SBarry Smith     }
7933ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
794f3fbd535SBarry Smith     for (i=1; i<n; i++) {
7953ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
7963ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
797f3fbd535SBarry Smith     }
798f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
799f3fbd535SBarry Smith       if (!mglevels[i]->b) {
800f3fbd535SBarry Smith         Vec *vec;
8012a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
802f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
8036bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
804f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
805f3fbd535SBarry Smith       }
806f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
807f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
808f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
8096bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
810f3fbd535SBarry Smith       }
811f3fbd535SBarry Smith       if (!mglevels[i]->x) {
812f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
813f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
8146bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
815f3fbd535SBarry Smith       }
816f3fbd535SBarry Smith     }
817f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
818f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
819f3fbd535SBarry Smith       Vec *vec;
8202a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
821f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
8226bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
823f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
824f3fbd535SBarry Smith     }
825f3fbd535SBarry Smith   }
826f3fbd535SBarry Smith 
82798e04984SBarry Smith   if (pc->dm) {
82898e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
82998e04984SBarry Smith     for (i=0; i<n-1; i++) {
83098e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
83198e04984SBarry Smith     }
83298e04984SBarry Smith   }
833f3fbd535SBarry Smith 
834f3fbd535SBarry Smith   for (i=1; i<n; i++) {
8352a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
836f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
837f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
838f3fbd535SBarry Smith     }
83963e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
840f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
841*c0decd05SBarry Smith     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
842899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
843899639b0SHong Zhang     }
84463e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
845d42688cbSBarry Smith     if (!mglevels[i]->residual) {
846d42688cbSBarry Smith       Mat mat;
84713842ffbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
84854b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
849d42688cbSBarry Smith     }
850f3fbd535SBarry Smith   }
851f3fbd535SBarry Smith   for (i=1; i<n; i++) {
852f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
853f3fbd535SBarry Smith       Mat downmat,downpmat;
854f3fbd535SBarry Smith 
855f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
8560298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
857f3fbd535SBarry Smith       if (!opsset) {
85823ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
85923ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
860f3fbd535SBarry Smith       }
861f3fbd535SBarry Smith 
862f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
86363e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
864f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
865*c0decd05SBarry Smith       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
866899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
867899639b0SHong Zhang       }
86863e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
869f3fbd535SBarry Smith     }
870f3fbd535SBarry Smith   }
871f3fbd535SBarry Smith 
87263e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
873f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
874*c0decd05SBarry Smith   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
875899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
876899639b0SHong Zhang   }
87763e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
878f3fbd535SBarry Smith 
879f3fbd535SBarry Smith   /*
880f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
881e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
882f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
883f3fbd535SBarry Smith 
884f3fbd535SBarry Smith    Only support one or the other at the same time.
885f3fbd535SBarry Smith   */
886f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
887c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
888ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
889f3fbd535SBarry Smith   dump = PETSC_FALSE;
890f3fbd535SBarry Smith #endif
891c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
892ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
893f3fbd535SBarry Smith 
894f3fbd535SBarry Smith   if (viewer) {
895f3fbd535SBarry Smith     for (i=1; i<n; i++) {
896f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
897f3fbd535SBarry Smith     }
898f3fbd535SBarry Smith     for (i=0; i<n; i++) {
899f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
900f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
901f3fbd535SBarry Smith     }
902f3fbd535SBarry Smith   }
903f3fbd535SBarry Smith   PetscFunctionReturn(0);
904f3fbd535SBarry Smith }
905f3fbd535SBarry Smith 
906f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
907f3fbd535SBarry Smith 
9084b9ad928SBarry Smith /*@
90997177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
9104b9ad928SBarry Smith 
9114b9ad928SBarry Smith    Not Collective
9124b9ad928SBarry Smith 
9134b9ad928SBarry Smith    Input Parameter:
9144b9ad928SBarry Smith .  pc - the preconditioner context
9154b9ad928SBarry Smith 
9164b9ad928SBarry Smith    Output parameter:
9174b9ad928SBarry Smith .  levels - the number of levels
9184b9ad928SBarry Smith 
9194b9ad928SBarry Smith    Level: advanced
9204b9ad928SBarry Smith 
9214b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
9224b9ad928SBarry Smith 
92397177400SBarry Smith .seealso: PCMGSetLevels()
9244b9ad928SBarry Smith @*/
9257087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
9264b9ad928SBarry Smith {
927f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9284b9ad928SBarry Smith 
9294b9ad928SBarry Smith   PetscFunctionBegin;
9300700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9314482741eSBarry Smith   PetscValidIntPointer(levels,2);
932f3fbd535SBarry Smith   *levels = mg->nlevels;
9334b9ad928SBarry Smith   PetscFunctionReturn(0);
9344b9ad928SBarry Smith }
9354b9ad928SBarry Smith 
9364b9ad928SBarry Smith /*@
93797177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
9384b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
9394b9ad928SBarry Smith 
940ad4df100SBarry Smith    Logically Collective on PC
9414b9ad928SBarry Smith 
9424b9ad928SBarry Smith    Input Parameters:
9434b9ad928SBarry Smith +  pc - the preconditioner context
9449dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
9459dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
9464b9ad928SBarry Smith 
9474b9ad928SBarry Smith    Options Database Key:
9484b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
9494b9ad928SBarry Smith    additive, full, kaskade
9504b9ad928SBarry Smith 
9514b9ad928SBarry Smith    Level: advanced
9524b9ad928SBarry Smith 
9534b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
9544b9ad928SBarry Smith 
95597177400SBarry Smith .seealso: PCMGSetLevels()
9564b9ad928SBarry Smith @*/
9577087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
9584b9ad928SBarry Smith {
959f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9604b9ad928SBarry Smith 
9614b9ad928SBarry Smith   PetscFunctionBegin;
9620700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
963c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
96431567311SBarry Smith   mg->am = form;
9659dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
966c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
967c60c7ad4SBarry Smith   PetscFunctionReturn(0);
968c60c7ad4SBarry Smith }
969c60c7ad4SBarry Smith 
970c60c7ad4SBarry Smith /*@
971c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
972c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
973c60c7ad4SBarry Smith 
974c60c7ad4SBarry Smith    Logically Collective on PC
975c60c7ad4SBarry Smith 
976c60c7ad4SBarry Smith    Input Parameter:
977c60c7ad4SBarry Smith .  pc - the preconditioner context
978c60c7ad4SBarry Smith 
979c60c7ad4SBarry Smith    Output Parameter:
980c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
981c60c7ad4SBarry Smith 
982c60c7ad4SBarry Smith 
983c60c7ad4SBarry Smith    Level: advanced
984c60c7ad4SBarry Smith 
985c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
986c60c7ad4SBarry Smith 
987c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
988c60c7ad4SBarry Smith @*/
989c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
990c60c7ad4SBarry Smith {
991c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
992c60c7ad4SBarry Smith 
993c60c7ad4SBarry Smith   PetscFunctionBegin;
994c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
995c60c7ad4SBarry Smith   *type = mg->am;
9964b9ad928SBarry Smith   PetscFunctionReturn(0);
9974b9ad928SBarry Smith }
9984b9ad928SBarry Smith 
9994b9ad928SBarry Smith /*@
10000d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
10014b9ad928SBarry Smith    complicated cycling.
10024b9ad928SBarry Smith 
1003ad4df100SBarry Smith    Logically Collective on PC
10044b9ad928SBarry Smith 
10054b9ad928SBarry Smith    Input Parameters:
1006c2be2410SBarry Smith +  pc - the multigrid context
1007c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
10084b9ad928SBarry Smith 
10094b9ad928SBarry Smith    Options Database Key:
1010c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
10114b9ad928SBarry Smith 
10124b9ad928SBarry Smith    Level: advanced
10134b9ad928SBarry Smith 
10144b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10154b9ad928SBarry Smith 
10160d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
10174b9ad928SBarry Smith @*/
10187087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
10194b9ad928SBarry Smith {
1020f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
1021f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
102279416396SBarry Smith   PetscInt     i,levels;
10234b9ad928SBarry Smith 
10244b9ad928SBarry Smith   PetscFunctionBegin;
10250700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
102669fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
10271a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1028f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10292fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
10304b9ad928SBarry Smith   PetscFunctionReturn(0);
10314b9ad928SBarry Smith }
10324b9ad928SBarry Smith 
10338cc2d5dfSBarry Smith /*@
10348cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
10358cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
10368cc2d5dfSBarry Smith 
1037ad4df100SBarry Smith    Logically Collective on PC
10388cc2d5dfSBarry Smith 
10398cc2d5dfSBarry Smith    Input Parameters:
10408cc2d5dfSBarry Smith +  pc - the multigrid context
10418cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
10428cc2d5dfSBarry Smith 
10438cc2d5dfSBarry Smith    Options Database Key:
1044e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
10458cc2d5dfSBarry Smith 
10468cc2d5dfSBarry Smith    Level: advanced
10478cc2d5dfSBarry Smith 
104895452b02SPatrick Sanan    Notes:
104995452b02SPatrick Sanan     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
10508cc2d5dfSBarry Smith 
10518cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10528cc2d5dfSBarry Smith 
10538cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
10548cc2d5dfSBarry Smith @*/
10557087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
10568cc2d5dfSBarry Smith {
1057f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
10588cc2d5dfSBarry Smith 
10598cc2d5dfSBarry Smith   PetscFunctionBegin;
10600700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1061c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
10622a21e185SBarry Smith   mg->cyclesperpcapply = n;
10638cc2d5dfSBarry Smith   PetscFunctionReturn(0);
10648cc2d5dfSBarry Smith }
10658cc2d5dfSBarry Smith 
10662134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1067fb15c04eSBarry Smith {
1068fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1069fb15c04eSBarry Smith 
1070fb15c04eSBarry Smith   PetscFunctionBegin;
10712134b1e4SBarry Smith   mg->galerkin = use;
1072fb15c04eSBarry Smith   PetscFunctionReturn(0);
1073fb15c04eSBarry Smith }
1074fb15c04eSBarry Smith 
1075c2be2410SBarry Smith /*@
107697177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
107782b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1078c2be2410SBarry Smith 
1079ad4df100SBarry Smith    Logically Collective on PC
1080c2be2410SBarry Smith 
1081c2be2410SBarry Smith    Input Parameters:
1082c91913d3SJed Brown +  pc - the multigrid context
10832134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1084c2be2410SBarry Smith 
1085c2be2410SBarry Smith    Options Database Key:
10862134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1087c2be2410SBarry Smith 
1088c2be2410SBarry Smith    Level: intermediate
1089c2be2410SBarry Smith 
109095452b02SPatrick Sanan    Notes:
109195452b02SPatrick Sanan     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
10921c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
10931c1aac46SBarry Smith 
1094c2be2410SBarry Smith .keywords: MG, set, Galerkin
1095c2be2410SBarry Smith 
10962134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
10973fc8bf9cSBarry Smith 
1098c2be2410SBarry Smith @*/
10992134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1100c2be2410SBarry Smith {
1101fb15c04eSBarry Smith   PetscErrorCode ierr;
1102c2be2410SBarry Smith 
1103c2be2410SBarry Smith   PetscFunctionBegin;
11040700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11052134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1106c2be2410SBarry Smith   PetscFunctionReturn(0);
1107c2be2410SBarry Smith }
1108c2be2410SBarry Smith 
11093fc8bf9cSBarry Smith /*@
11103fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
111182b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
11123fc8bf9cSBarry Smith 
11133fc8bf9cSBarry Smith    Not Collective
11143fc8bf9cSBarry Smith 
11153fc8bf9cSBarry Smith    Input Parameter:
11163fc8bf9cSBarry Smith .  pc - the multigrid context
11173fc8bf9cSBarry Smith 
11183fc8bf9cSBarry Smith    Output Parameter:
11192134b1e4SBarry Smith .  galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL
11203fc8bf9cSBarry Smith 
11213fc8bf9cSBarry Smith    Level: intermediate
11223fc8bf9cSBarry Smith 
11233fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
11243fc8bf9cSBarry Smith 
11252134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
11263fc8bf9cSBarry Smith 
11273fc8bf9cSBarry Smith @*/
11282134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
11293fc8bf9cSBarry Smith {
1130f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
11313fc8bf9cSBarry Smith 
11323fc8bf9cSBarry Smith   PetscFunctionBegin;
11330700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11342134b1e4SBarry Smith   *galerkin = mg->galerkin;
11353fc8bf9cSBarry Smith   PetscFunctionReturn(0);
11363fc8bf9cSBarry Smith }
11373fc8bf9cSBarry Smith 
11384b9ad928SBarry Smith /*@
113906792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1140710315b6SLawrence Mitchell    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1141710315b6SLawrence Mitchell    pre- and post-smoothing steps.
114206792cafSBarry Smith 
114306792cafSBarry Smith    Logically Collective on PC
114406792cafSBarry Smith 
114506792cafSBarry Smith    Input Parameters:
114606792cafSBarry Smith +  mg - the multigrid context
114706792cafSBarry Smith -  n - the number of smoothing steps
114806792cafSBarry Smith 
114906792cafSBarry Smith    Options Database Key:
115006792cafSBarry Smith +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
115106792cafSBarry Smith 
115206792cafSBarry Smith    Level: advanced
115306792cafSBarry Smith 
115495452b02SPatrick Sanan    Notes:
115595452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
115606792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
115706792cafSBarry Smith 
115806792cafSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
115906792cafSBarry Smith 
1160710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp()
116106792cafSBarry Smith @*/
116206792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
116306792cafSBarry Smith {
116406792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
116506792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
116606792cafSBarry Smith   PetscErrorCode ierr;
116706792cafSBarry Smith   PetscInt       i,levels;
116806792cafSBarry Smith 
116906792cafSBarry Smith   PetscFunctionBegin;
117006792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
117106792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
11721a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
117306792cafSBarry Smith   levels = mglevels[0]->levels;
117406792cafSBarry Smith 
117506792cafSBarry Smith   for (i=1; i<levels; i++) {
117606792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
117706792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
117806792cafSBarry Smith     mg->default_smoothu = n;
117906792cafSBarry Smith     mg->default_smoothd = n;
118006792cafSBarry Smith   }
118106792cafSBarry Smith   PetscFunctionReturn(0);
118206792cafSBarry Smith }
118306792cafSBarry Smith 
1184f442ab6aSBarry Smith /*@
1185f442ab6aSBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels
1186710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1187f442ab6aSBarry Smith 
1188f442ab6aSBarry Smith    Logically Collective on PC
1189f442ab6aSBarry Smith 
1190f442ab6aSBarry Smith    Input Parameters:
1191f442ab6aSBarry Smith .  pc - the preconditioner context
1192f442ab6aSBarry Smith 
1193f442ab6aSBarry Smith    Options Database Key:
1194f442ab6aSBarry Smith .  -pc_mg_distinct_smoothup
1195f442ab6aSBarry Smith 
1196f442ab6aSBarry Smith    Level: advanced
1197f442ab6aSBarry Smith 
119895452b02SPatrick Sanan    Notes:
119995452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
1200f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1201f442ab6aSBarry Smith 
1202f442ab6aSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1203f442ab6aSBarry Smith 
1204710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth()
1205f442ab6aSBarry Smith @*/
1206f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1207f442ab6aSBarry Smith {
1208f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1209f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1210f442ab6aSBarry Smith   PetscErrorCode ierr;
1211f442ab6aSBarry Smith   PetscInt       i,levels;
1212f442ab6aSBarry Smith   KSP            subksp;
1213f442ab6aSBarry Smith 
1214f442ab6aSBarry Smith   PetscFunctionBegin;
1215f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1216f442ab6aSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1217f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1218f442ab6aSBarry Smith 
1219f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1220710315b6SLawrence Mitchell     const char *prefix = NULL;
1221f442ab6aSBarry Smith     /* make sure smoother up and down are different */
1222f442ab6aSBarry Smith     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1223710315b6SLawrence Mitchell     ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr);
1224710315b6SLawrence Mitchell     ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr);
1225f442ab6aSBarry Smith     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1226f442ab6aSBarry Smith   }
1227f442ab6aSBarry Smith   PetscFunctionReturn(0);
1228f442ab6aSBarry Smith }
1229f442ab6aSBarry Smith 
12304b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
12314b9ad928SBarry Smith 
12323b09bd56SBarry Smith /*MC
1233ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
12343b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
12353b09bd56SBarry Smith 
12363b09bd56SBarry Smith    Options Database Keys:
12373b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1238391689abSStefano Zampini .  -pc_mg_cycle_type <v,w> - provide the cycle desired
12398c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
12403b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1241710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
12422134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
12438cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
12448cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1245e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
12468cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
12478cf6037eSBarry Smith                         to the binary output file called binaryoutput
12483b09bd56SBarry Smith 
124995452b02SPatrick Sanan    Notes:
12500b3c753eSRichard Tran Mills     If one uses a Krylov method such GMRES or CG as the smoother then one must use KSPFGMRES, KSPGCR, or KSPRICHARDSON as the outer Krylov method
12513b09bd56SBarry Smith 
12528cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
12538cf6037eSBarry Smith 
125423067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
125523067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
125623067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
125723067569SBarry Smith        residual is computed at the end of each cycle.
125823067569SBarry Smith 
12593b09bd56SBarry Smith    Level: intermediate
12603b09bd56SBarry Smith 
12618f87f92bSBarry Smith    Concepts: multigrid/multilevel
12623b09bd56SBarry Smith 
12638cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1264710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1265710315b6SLawrence Mitchell            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
126697177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12670d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12683b09bd56SBarry Smith M*/
12693b09bd56SBarry Smith 
12708cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
12714b9ad928SBarry Smith {
1272f3fbd535SBarry Smith   PC_MG          *mg;
1273f3fbd535SBarry Smith   PetscErrorCode ierr;
1274f3fbd535SBarry Smith 
12754b9ad928SBarry Smith   PetscFunctionBegin;
1276b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1277f3fbd535SBarry Smith   pc->data     = (void*)mg;
1278f3fbd535SBarry Smith   mg->nlevels  = -1;
127910eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
12802134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1281f3fbd535SBarry Smith 
128237a44384SMark Adams   pc->useAmat = PETSC_TRUE;
128337a44384SMark Adams 
12844b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
12854b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1286a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
12874b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
12884b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
12894b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1290fb15c04eSBarry Smith 
1291fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
12924b9ad928SBarry Smith   PetscFunctionReturn(0);
12934b9ad928SBarry Smith }
1294