xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 41fce8fd74210c22031f6993d13f4ee58a42f2d0)
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>
74b9ad928SBarry Smith 
831567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
94b9ad928SBarry Smith {
1031567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1131567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
126849ba73SBarry Smith   PetscErrorCode ierr;
1331567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
14a0808db4SHong Zhang   PC             subpc;
15a0808db4SHong Zhang   PCFailedReason pcreason;
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 */
20a0808db4SHong Zhang   ierr = KSPGetPC(mglevels->smoothd,&subpc);CHKERRQ(ierr);
21a0808db4SHong Zhang   ierr = PCGetSetUpFailedReason(subpc,&pcreason);CHKERRQ(ierr);
22a0808db4SHong Zhang   if (pcreason) {
23a0808db4SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
24a0808db4SHong Zhang   }
2563e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2631567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2763e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2831567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2963e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
304b9ad928SBarry Smith 
314b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
3231567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
334b9ad928SBarry Smith       PetscReal rnorm;
3431567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
354b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3670441072SBarry Smith         if (rnorm < mg->abstol) {
374d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
3857622a8eSBarry 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);
394b9ad928SBarry Smith         } else {
404d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
4157622a8eSBarry 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);
424b9ad928SBarry Smith         }
434b9ad928SBarry Smith         PetscFunctionReturn(0);
444b9ad928SBarry Smith       }
454b9ad928SBarry Smith     }
464b9ad928SBarry Smith 
4731567311SBarry Smith     mgc = *(mglevelsin - 1);
4863e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4931567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
5063e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
51efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
524b9ad928SBarry Smith     while (cycles--) {
5331567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
544b9ad928SBarry Smith     }
5563e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5631567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5763e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5863e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5931567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
6063e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
614b9ad928SBarry Smith   }
624b9ad928SBarry Smith   PetscFunctionReturn(0);
634b9ad928SBarry Smith }
644b9ad928SBarry Smith 
65ace3abfcSBarry 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)
664b9ad928SBarry Smith {
67f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
68f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
69dfbe8321SBarry Smith   PetscErrorCode ierr;
70f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
714b9ad928SBarry Smith 
724b9ad928SBarry Smith   PetscFunctionBegin;
73a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
74a762d673SBarry Smith   for (i=0; i<levels; i++) {
75a762d673SBarry Smith     if (!mglevels[i]->A) {
76a762d673SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
77a762d673SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
78a762d673SBarry Smith     }
79a762d673SBarry Smith   }
80f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
81f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
824b9ad928SBarry Smith 
8331567311SBarry Smith   mg->rtol   = rtol;
8431567311SBarry Smith   mg->abstol = abstol;
8531567311SBarry Smith   mg->dtol   = dtol;
864b9ad928SBarry Smith   if (rtol) {
874b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
884b9ad928SBarry Smith     PetscReal rnorm;
897319c654SBarry Smith     if (zeroguess) {
907319c654SBarry Smith       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
917319c654SBarry Smith     } else {
92f3fbd535SBarry Smith       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
934b9ad928SBarry Smith       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
947319c654SBarry Smith     }
9531567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
962fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
972fa5cd67SKarl Rupp   else mg->ttol = 0.0;
984b9ad928SBarry Smith 
994d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
100336babb1SJed Brown      stop prematurely due to small residual */
1014d0a8057SBarry Smith   for (i=1; i<levels; i++) {
102f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
103f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
10423067569SBarry Smith       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
10523067569SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
106f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1074b9ad928SBarry Smith     }
1084d0a8057SBarry Smith   }
1094d0a8057SBarry Smith 
1104d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1114d0a8057SBarry Smith   for (i=0; i<its; i++) {
112f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1134d0a8057SBarry Smith     if (*reason) break;
1144d0a8057SBarry Smith   }
1154d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1164d0a8057SBarry Smith   *outits = i;
1174b9ad928SBarry Smith   PetscFunctionReturn(0);
1184b9ad928SBarry Smith }
1194b9ad928SBarry Smith 
1203aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1213aeaf226SBarry Smith {
1223aeaf226SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1233aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1243aeaf226SBarry Smith   PetscErrorCode ierr;
1253aeaf226SBarry Smith   PetscInt       i,n;
1263aeaf226SBarry Smith 
1273aeaf226SBarry Smith   PetscFunctionBegin;
1283aeaf226SBarry Smith   if (mglevels) {
1293aeaf226SBarry Smith     n = mglevels[0]->levels;
1303aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1313aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1323aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1333aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
1343aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1353aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
13673250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1373aeaf226SBarry Smith     }
1383aeaf226SBarry Smith 
1393aeaf226SBarry Smith     for (i=0; i<n; i++) {
1403aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1413aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1423aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1433aeaf226SBarry Smith       }
1443aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1453aeaf226SBarry Smith     }
1463aeaf226SBarry Smith   }
1473aeaf226SBarry Smith   PetscFunctionReturn(0);
1483aeaf226SBarry Smith }
1493aeaf226SBarry Smith 
1504b9ad928SBarry Smith /*@C
15197177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1524b9ad928SBarry Smith    Must be called before any other MG routine.
1534b9ad928SBarry Smith 
154ad4df100SBarry Smith    Logically Collective on PC
1554b9ad928SBarry Smith 
1564b9ad928SBarry Smith    Input Parameters:
1574b9ad928SBarry Smith +  pc - the preconditioner context
1584b9ad928SBarry Smith .  levels - the number of levels
1594b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1602bf68e3eSBarry Smith            on smaller sets of processors.
1614b9ad928SBarry Smith 
1624b9ad928SBarry Smith    Level: intermediate
1634b9ad928SBarry Smith 
1644b9ad928SBarry Smith    Notes:
1654b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1664b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1674b9ad928SBarry Smith 
1684b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1694b9ad928SBarry Smith 
17097177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1714b9ad928SBarry Smith @*/
1727087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1734b9ad928SBarry Smith {
174dfbe8321SBarry Smith   PetscErrorCode ierr;
175f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
176ce94432eSBarry Smith   MPI_Comm       comm;
1773aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
17810eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
17910167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
180f3fbd535SBarry Smith   PetscInt       i;
181f3fbd535SBarry Smith   PetscMPIInt    size;
182f3fbd535SBarry Smith   const char     *prefix;
183f3fbd535SBarry Smith   PC             ipc;
1843aeaf226SBarry Smith   PetscInt       n;
1854b9ad928SBarry Smith 
1864b9ad928SBarry Smith   PetscFunctionBegin;
1870700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
188548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
189ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
190548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1913aeaf226SBarry Smith   if (mglevels) {
19210eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
1933aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
1943aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
1953aeaf226SBarry Smith     n    = mglevels[0]->levels;
1963aeaf226SBarry Smith     for (i=0; i<n; i++) {
1973aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1983aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
1993aeaf226SBarry Smith       }
2003aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
2013aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
2023aeaf226SBarry Smith     }
2033aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
2043aeaf226SBarry Smith   }
205f3fbd535SBarry Smith 
206f3fbd535SBarry Smith   mg->nlevels = levels;
207f3fbd535SBarry Smith 
208785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
2093bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
210f3fbd535SBarry Smith 
211f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
212f3fbd535SBarry Smith 
213a9db87a2SMatthew G Knepley   mg->stageApply = 0;
214f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
215b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
2162fa5cd67SKarl Rupp 
217f3fbd535SBarry Smith     mglevels[i]->level               = i;
218f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
21910eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
220336babb1SJed Brown     mg->default_smoothu              = 2;
221336babb1SJed Brown     mg->default_smoothd              = 2;
22263e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
22363e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
22463e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
22563e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
226f3fbd535SBarry Smith 
227f3fbd535SBarry Smith     if (comms) comm = comms[i];
228f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
229422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
2300ee9a668SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
2310ee9a668SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
2320ee9a668SBarry Smith     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
2330ee9a668SBarry Smith     if (i || levels == 1) {
2340ee9a668SBarry Smith       char tprefix[128];
2350ee9a668SBarry Smith 
236336babb1SJed Brown       ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2370059c7bdSJed Brown       ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
238336babb1SJed Brown       ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
239336babb1SJed Brown       ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
240336babb1SJed Brown       ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
2410ee9a668SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
242f3fbd535SBarry Smith 
2430ee9a668SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
2440ee9a668SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
2450ee9a668SBarry Smith     } else {
246f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
247f3fbd535SBarry Smith 
2489dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
249f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
250f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
251f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
252f3fbd535SBarry Smith       if (size > 1) {
253f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
254f3fbd535SBarry Smith       } else {
255f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
256f3fbd535SBarry Smith       }
257753b7fb9SBarry Smith       ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
258f3fbd535SBarry Smith     }
2593bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
2602fa5cd67SKarl Rupp 
261f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
26231567311SBarry Smith     mg->rtol             = 0.0;
26331567311SBarry Smith     mg->abstol           = 0.0;
26431567311SBarry Smith     mg->dtol             = 0.0;
26531567311SBarry Smith     mg->ttol             = 0.0;
26631567311SBarry Smith     mg->cyclesperpcapply = 1;
267f3fbd535SBarry Smith   }
268f3fbd535SBarry Smith   mg->levels = mglevels;
26910eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2704b9ad928SBarry Smith   PetscFunctionReturn(0);
2714b9ad928SBarry Smith }
2724b9ad928SBarry Smith 
273c07bf074SBarry Smith 
274c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
275c07bf074SBarry Smith {
276c07bf074SBarry Smith   PetscErrorCode ierr;
277a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
278a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
279a06653b4SBarry Smith   PetscInt       i,n;
280c07bf074SBarry Smith 
281c07bf074SBarry Smith   PetscFunctionBegin;
282a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
283a06653b4SBarry Smith   if (mglevels) {
284a06653b4SBarry Smith     n = mglevels[0]->levels;
285a06653b4SBarry Smith     for (i=0; i<n; i++) {
286a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2876bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
288a06653b4SBarry Smith       }
2896bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
290a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
291a06653b4SBarry Smith     }
292a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
293a06653b4SBarry Smith   }
294c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
295f3fbd535SBarry Smith   PetscFunctionReturn(0);
296f3fbd535SBarry Smith }
297f3fbd535SBarry Smith 
298f3fbd535SBarry Smith 
299f3fbd535SBarry Smith 
30009573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
30109573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
30209573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
303f3fbd535SBarry Smith 
304f3fbd535SBarry Smith /*
305f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
306f3fbd535SBarry Smith              or full cycle of multigrid.
307f3fbd535SBarry Smith 
308f3fbd535SBarry Smith   Note:
309f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
310f3fbd535SBarry Smith */
311f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
312f3fbd535SBarry Smith {
313f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
314f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
315f3fbd535SBarry Smith   PetscErrorCode ierr;
316f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
317f3fbd535SBarry Smith 
318f3fbd535SBarry Smith   PetscFunctionBegin;
319a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
320e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
321298cc208SBarry Smith   for (i=0; i<levels; i++) {
322e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
32323ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
324298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
325e1d8e5deSBarry Smith     }
326e1d8e5deSBarry Smith   }
327e1d8e5deSBarry Smith 
328f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
329f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
33031567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
331f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
33231567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
3330298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
334f3fbd535SBarry Smith     }
3352fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
33631567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
3372fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
33831567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
3392fa5cd67SKarl Rupp   } 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 
3474416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
348f3fbd535SBarry Smith {
349f3fbd535SBarry Smith   PetscErrorCode   ierr;
350f3fbd535SBarry Smith   PetscInt         m,levels = 1,cycles;
3512134b1e4SBarry Smith   PetscBool        flg;
352f3fbd535SBarry Smith   PC_MG            *mg        = (PC_MG*)pc->data;
3533d3eaba7SBarry Smith   PC_MG_Levels     **mglevels;
354f3fbd535SBarry Smith   PCMGType         mgtype;
355f3fbd535SBarry Smith   PCMGCycleType    mgctype;
3562134b1e4SBarry Smith   PCMGGalerkinType gtype;
357f3fbd535SBarry Smith 
358f3fbd535SBarry Smith   PetscFunctionBegin;
359e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
3603aeaf226SBarry Smith   if (!mg->levels) {
361f3fbd535SBarry Smith     ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
362eb3f98d2SBarry Smith     if (!flg && pc->dm) {
363eb3f98d2SBarry Smith       ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
364eb3f98d2SBarry Smith       levels++;
3653aeaf226SBarry Smith       mg->usedmfornumberoflevels = PETSC_TRUE;
366eb3f98d2SBarry Smith     }
3670298fd71SBarry Smith     ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
368f3fbd535SBarry Smith   }
3693aeaf226SBarry Smith   mglevels = mg->levels;
3703aeaf226SBarry Smith 
371f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
372f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
373f3fbd535SBarry Smith   if (flg) {
374f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
3752fa5cd67SKarl Rupp   }
3762134b1e4SBarry Smith   gtype = mg->galerkin;
3772134b1e4SBarry Smith   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
3782134b1e4SBarry Smith   if (flg) {
3792134b1e4SBarry Smith     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
380f3fbd535SBarry Smith   }
38194ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);CHKERRQ(ierr);
382f3fbd535SBarry Smith   if (flg) {
383f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
384f3fbd535SBarry Smith   }
38594ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);CHKERRQ(ierr);
386f3fbd535SBarry Smith   if (flg) {
387f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
388f3fbd535SBarry Smith   }
38931567311SBarry Smith   mgtype = mg->am;
390f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
391f3fbd535SBarry Smith   if (flg) {
392f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
393f3fbd535SBarry Smith   }
39431567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
3955363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
396f3fbd535SBarry Smith     if (flg) {
397f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
398f3fbd535SBarry Smith     }
399f3fbd535SBarry Smith   }
400f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4010298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
402f3fbd535SBarry Smith   if (flg) {
403f3fbd535SBarry Smith     PetscInt i;
404f3fbd535SBarry Smith     char     eventname[128];
405ce94432eSBarry Smith     if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
406f3fbd535SBarry Smith     levels = mglevels[0]->levels;
407f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
408f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
40963e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
410f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
41163e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
412f3fbd535SBarry Smith       if (i) {
413f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
41463e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
415f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
41663e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
417f3fbd535SBarry Smith       }
418f3fbd535SBarry Smith     }
419eec35531SMatthew G Knepley 
420ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
421eec35531SMatthew G Knepley     {
422eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
423eec35531SMatthew G Knepley       PetscStageLog stageLog;
424eec35531SMatthew G Knepley       PetscInt      st;
425eec35531SMatthew G Knepley 
426eec35531SMatthew G Knepley       PetscFunctionBegin;
427eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
428eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
429eec35531SMatthew G Knepley         PetscBool same;
430eec35531SMatthew G Knepley 
431eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4322fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
433eec35531SMatthew G Knepley       }
434eec35531SMatthew G Knepley       if (!mg->stageApply) {
435eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
436eec35531SMatthew G Knepley       }
437eec35531SMatthew G Knepley     }
438ec60ca73SSatish Balay #endif
439f3fbd535SBarry Smith   }
440f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
441f3fbd535SBarry Smith   PetscFunctionReturn(0);
442f3fbd535SBarry Smith }
443f3fbd535SBarry Smith 
4446a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4456a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
4462134b1e4SBarry Smith const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0};
447f3fbd535SBarry Smith 
4489804daf3SBarry Smith #include <petscdraw.h>
449f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
450f3fbd535SBarry Smith {
451f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
452f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
453f3fbd535SBarry Smith   PetscErrorCode ierr;
454e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
455219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
456f3fbd535SBarry Smith 
457f3fbd535SBarry Smith   PetscFunctionBegin;
458251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4595b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
460219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
461f3fbd535SBarry Smith   if (iascii) {
462e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
463efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);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     }
4672134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
468f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4692134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
4702134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
4712134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
4722134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
4732134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
4742134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
4754f66f45eSBarry Smith     } else {
4764f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
477f3fbd535SBarry Smith     }
4785adeb434SBarry Smith     if (mg->view){
4795adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
4805adeb434SBarry Smith     }
481f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
482f3fbd535SBarry Smith       if (!i) {
48389cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
484f3fbd535SBarry Smith       } else {
48589cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
486f3fbd535SBarry Smith       }
487f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
488f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
489f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
490f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
491f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
492f3fbd535SBarry Smith       } else if (i) {
49389cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
494f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
495f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
496f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
497f3fbd535SBarry Smith       }
498f3fbd535SBarry Smith     }
4995b0b0462SBarry Smith   } else if (isbinary) {
5005b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5015b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5025b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5035b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5045b0b0462SBarry Smith       }
5055b0b0462SBarry Smith     }
506219b1045SBarry Smith   } else if (isdraw) {
507219b1045SBarry Smith     PetscDraw draw;
508b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
509219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
510219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5110298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
512219b1045SBarry Smith     bottom = y - th;
513219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
514b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
515219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
516219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
517219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
518b4375e8dSPeter Brune       } else {
519b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
520b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
521b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
522b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
523b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
524b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
525b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
526b4375e8dSPeter Brune       }
5270298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5281cd381d6SBarry Smith       bottom -= th;
529219b1045SBarry Smith     }
5305b0b0462SBarry Smith   }
531f3fbd535SBarry Smith   PetscFunctionReturn(0);
532f3fbd535SBarry Smith }
533f3fbd535SBarry Smith 
534af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
535af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
536cab2e9ccSBarry Smith 
537f3fbd535SBarry Smith /*
538f3fbd535SBarry Smith     Calls setup for the KSP on each level
539f3fbd535SBarry Smith */
540f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
541f3fbd535SBarry Smith {
542f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
543f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
544f3fbd535SBarry Smith   PetscErrorCode ierr;
545f3fbd535SBarry Smith   PetscInt       i,n = mglevels[0]->levels;
54698e04984SBarry Smith   PC             cpc;
5473db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
548f3fbd535SBarry Smith   Mat            dA,dB;
549f3fbd535SBarry Smith   Vec            tvec;
550218a07d4SBarry Smith   DM             *dms;
551649052a6SBarry Smith   PetscViewer    viewer = 0;
5522134b1e4SBarry Smith   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
553f3fbd535SBarry Smith 
554f3fbd535SBarry Smith   PetscFunctionBegin;
55501bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5563aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5573aeaf226SBarry Smith     PetscInt levels;
5583aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5593aeaf226SBarry Smith     levels++;
5603aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5610298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5623aeaf226SBarry Smith       n        = levels;
5633aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5643aeaf226SBarry Smith       mglevels =  mg->levels;
5653aeaf226SBarry Smith     }
5663aeaf226SBarry Smith   }
56798e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5683aeaf226SBarry Smith 
569f3fbd535SBarry Smith 
570f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
571f3fbd535SBarry Smith   /* so use those from global PC */
572f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
5730298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
57498e04984SBarry Smith   if (opsset) {
57598e04984SBarry Smith     Mat mmat;
57623ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
57798e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
57898e04984SBarry Smith   }
5794a5f07a5SMark F. Adams 
5804a5f07a5SMark F. Adams   if (!opsset) {
58171b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
5824a5f07a5SMark F. Adams     if(use_amat){
583f3fbd535SBarry 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);
58423ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
585f3fbd535SBarry Smith     }
5864a5f07a5SMark F. Adams     else {
5874a5f07a5SMark 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);
58823ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
5894a5f07a5SMark F. Adams     }
5904a5f07a5SMark F. Adams   }
591f3fbd535SBarry Smith 
5929c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
5939c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
5949c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
5959c7ffb39SBarry Smith       continue;
5969c7ffb39SBarry Smith     }
5979c7ffb39SBarry Smith   }
5982134b1e4SBarry Smith 
5992134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
6002134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
6012134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
6022134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
6032134b1e4SBarry Smith   }
6042134b1e4SBarry Smith 
6052134b1e4SBarry Smith 
6069c7ffb39SBarry Smith   /*
6079c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6089c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6099c7ffb39SBarry Smith   */
6102134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
6112d2b81a6SBarry Smith 	/* construct the interpolation from the DMs */
6122e499ae9SBarry Smith     Mat p;
61373250ac0SBarry Smith     Vec rscale;
614785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
615218a07d4SBarry Smith     dms[n-1] = pc->dm;
616ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
617ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
618*41fce8fdSHong Zhang 	/*
619*41fce8fdSHong Zhang 	   Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
620*41fce8fdSHong Zhang 	   Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
621*41fce8fdSHong Zhang 	   But it is safe to use -dm_mat_type.
622*41fce8fdSHong Zhang 
623*41fce8fdSHong Zhang 	   The mat type should not be hardcoded like this, we need to find a better way.
624*41fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
625*41fce8fdSHong Zhang     */
626218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
627942e3340SBarry Smith       DMKSP     kdm;
6283ad4599aSBarry Smith       PetscBool dmhasrestrict;
6293c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6302134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
631942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
632d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
633d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6340298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6350298fd71SBarry Smith       kdm->rhsctx          = NULL;
63624c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
637e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6382d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
63900ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
64073250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6416bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
642218a07d4SBarry Smith       }
6433ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
6443ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
6453ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
6463ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
6473ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
6483ad4599aSBarry Smith       }
64924c3aa18SBarry Smith     }
6502d2b81a6SBarry Smith 
651ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
6522d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
653ce4cda84SJed Brown   }
654cab2e9ccSBarry Smith 
655ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
6562134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
657cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
658cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
659218a07d4SBarry Smith   }
660218a07d4SBarry Smith 
6612134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
6622134b1e4SBarry Smith     Mat       A,B;
6632134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
6642134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
6652134b1e4SBarry Smith 
6662134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
6672134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
6682134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
669f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
670b9367914SBarry 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");
671b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
672b9367914SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
673b9367914SBarry Smith       }
674b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
675b9367914SBarry Smith         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
676b9367914SBarry Smith       }
6772134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
6782134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
6792134b1e4SBarry Smith       }
6802134b1e4SBarry Smith       if (doA) {
6812df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
6822134b1e4SBarry Smith       }
6832134b1e4SBarry Smith       if (doB) {
6842df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
685b9367914SBarry Smith       }
6862134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
6872134b1e4SBarry Smith       if (!doA && dAeqdB) {
6882134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
6892134b1e4SBarry Smith         A = B;
6902134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX ) {
6912134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
6922134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
693b9367914SBarry Smith       }
6942134b1e4SBarry Smith       if (!doB && dAeqdB) {
6952134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
6962134b1e4SBarry Smith         B = A;
6972134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
69823ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
6992134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
7007d537d92SJed Brown       }
7012134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
7022134b1e4SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
7032134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
7042134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
7052134b1e4SBarry Smith       }
7062134b1e4SBarry Smith       dA = A;
707f3fbd535SBarry Smith       dB = B;
708f3fbd535SBarry Smith     }
709f3fbd535SBarry Smith   }
7102134b1e4SBarry Smith   if (needRestricts && pc->dm && pc->dm->x) {
711cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
712cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
713c88c5224SJed Brown       Mat R;
714c88c5224SJed Brown       Vec rscale;
715cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
716cab2e9ccSBarry Smith         Vec *vecs;
7172a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
718cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
719cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
720cab2e9ccSBarry Smith       }
721c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
722c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
723c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
724c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
725cab2e9ccSBarry Smith     }
726f3fbd535SBarry Smith   }
7272134b1e4SBarry Smith   if (needRestricts && pc->dm) {
728caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
729ccceb50cSJed Brown       DM  dmfine,dmcoarse;
730caa4e7f2SJed Brown       Mat Restrict,Inject;
731caa4e7f2SJed Brown       Vec rscale;
732ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
733ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
734caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
735caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
7360298fd71SBarry Smith       Inject = NULL;      /* Callback should create it if it needs Injection */
737caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
738caa4e7f2SJed Brown     }
739caa4e7f2SJed Brown   }
740f3fbd535SBarry Smith 
741f3fbd535SBarry Smith   if (!pc->setupcalled) {
742f3fbd535SBarry Smith     for (i=0; i<n; i++) {
743f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
744f3fbd535SBarry Smith     }
745f3fbd535SBarry Smith     for (i=1; i<n; i++) {
746f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
747f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
748f3fbd535SBarry Smith       }
749f3fbd535SBarry Smith     }
7503ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
751f3fbd535SBarry Smith     for (i=1; i<n; i++) {
7523ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
7533ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
754f3fbd535SBarry Smith     }
755f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
756f3fbd535SBarry Smith       if (!mglevels[i]->b) {
757f3fbd535SBarry Smith         Vec *vec;
7582a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
759f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
7606bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
761f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
762f3fbd535SBarry Smith       }
763f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
764f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
765f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
7666bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
767f3fbd535SBarry Smith       }
768f3fbd535SBarry Smith       if (!mglevels[i]->x) {
769f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
770f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
7716bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
772f3fbd535SBarry Smith       }
773f3fbd535SBarry Smith     }
774f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
775f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
776f3fbd535SBarry Smith       Vec *vec;
7772a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
778f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
7796bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
780f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
781f3fbd535SBarry Smith     }
782f3fbd535SBarry Smith   }
783f3fbd535SBarry Smith 
78498e04984SBarry Smith   if (pc->dm) {
78598e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
78698e04984SBarry Smith     for (i=0; i<n-1; i++) {
78798e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
78898e04984SBarry Smith     }
78998e04984SBarry Smith   }
790f3fbd535SBarry Smith 
791f3fbd535SBarry Smith   for (i=1; i<n; i++) {
7922a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
793f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
794f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
795f3fbd535SBarry Smith     }
79663e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
797f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
798899639b0SHong Zhang     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
799899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
800899639b0SHong Zhang     }
80163e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
802d42688cbSBarry Smith     if (!mglevels[i]->residual) {
803d42688cbSBarry Smith       Mat mat;
80413842ffbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
80554b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
806d42688cbSBarry Smith     }
807f3fbd535SBarry Smith   }
808f3fbd535SBarry Smith   for (i=1; i<n; i++) {
809f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
810f3fbd535SBarry Smith       Mat downmat,downpmat;
811f3fbd535SBarry Smith 
812f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
8130298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
814f3fbd535SBarry Smith       if (!opsset) {
81523ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
81623ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
817f3fbd535SBarry Smith       }
818f3fbd535SBarry Smith 
819f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
82063e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
821f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
822899639b0SHong Zhang       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) {
823899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
824899639b0SHong Zhang       }
82563e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
826f3fbd535SBarry Smith     }
827f3fbd535SBarry Smith   }
828f3fbd535SBarry Smith 
82963e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
830f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
831899639b0SHong Zhang   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
832899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
833899639b0SHong Zhang   }
83463e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
835f3fbd535SBarry Smith 
836f3fbd535SBarry Smith   /*
837f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
838e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
839f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
840f3fbd535SBarry Smith 
841f3fbd535SBarry Smith    Only support one or the other at the same time.
842f3fbd535SBarry Smith   */
843f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
844c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
845ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
846f3fbd535SBarry Smith   dump = PETSC_FALSE;
847f3fbd535SBarry Smith #endif
848c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
849ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
850f3fbd535SBarry Smith 
851f3fbd535SBarry Smith   if (viewer) {
852f3fbd535SBarry Smith     for (i=1; i<n; i++) {
853f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
854f3fbd535SBarry Smith     }
855f3fbd535SBarry Smith     for (i=0; i<n; i++) {
856f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
857f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
858f3fbd535SBarry Smith     }
859f3fbd535SBarry Smith   }
860f3fbd535SBarry Smith   PetscFunctionReturn(0);
861f3fbd535SBarry Smith }
862f3fbd535SBarry Smith 
863f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
864f3fbd535SBarry Smith 
8654b9ad928SBarry Smith /*@
86697177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
8674b9ad928SBarry Smith 
8684b9ad928SBarry Smith    Not Collective
8694b9ad928SBarry Smith 
8704b9ad928SBarry Smith    Input Parameter:
8714b9ad928SBarry Smith .  pc - the preconditioner context
8724b9ad928SBarry Smith 
8734b9ad928SBarry Smith    Output parameter:
8744b9ad928SBarry Smith .  levels - the number of levels
8754b9ad928SBarry Smith 
8764b9ad928SBarry Smith    Level: advanced
8774b9ad928SBarry Smith 
8784b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8794b9ad928SBarry Smith 
88097177400SBarry Smith .seealso: PCMGSetLevels()
8814b9ad928SBarry Smith @*/
8827087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8834b9ad928SBarry Smith {
884f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8854b9ad928SBarry Smith 
8864b9ad928SBarry Smith   PetscFunctionBegin;
8870700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8884482741eSBarry Smith   PetscValidIntPointer(levels,2);
889f3fbd535SBarry Smith   *levels = mg->nlevels;
8904b9ad928SBarry Smith   PetscFunctionReturn(0);
8914b9ad928SBarry Smith }
8924b9ad928SBarry Smith 
8934b9ad928SBarry Smith /*@
89497177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8954b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8964b9ad928SBarry Smith 
897ad4df100SBarry Smith    Logically Collective on PC
8984b9ad928SBarry Smith 
8994b9ad928SBarry Smith    Input Parameters:
9004b9ad928SBarry Smith +  pc - the preconditioner context
9019dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
9029dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
9034b9ad928SBarry Smith 
9044b9ad928SBarry Smith    Options Database Key:
9054b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
9064b9ad928SBarry Smith    additive, full, kaskade
9074b9ad928SBarry Smith 
9084b9ad928SBarry Smith    Level: advanced
9094b9ad928SBarry Smith 
9104b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
9114b9ad928SBarry Smith 
91297177400SBarry Smith .seealso: PCMGSetLevels()
9134b9ad928SBarry Smith @*/
9147087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
9154b9ad928SBarry Smith {
916f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9174b9ad928SBarry Smith 
9184b9ad928SBarry Smith   PetscFunctionBegin;
9190700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
920c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
92131567311SBarry Smith   mg->am = form;
9229dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
923c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
924c60c7ad4SBarry Smith   PetscFunctionReturn(0);
925c60c7ad4SBarry Smith }
926c60c7ad4SBarry Smith 
927c60c7ad4SBarry Smith /*@
928c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
929c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
930c60c7ad4SBarry Smith 
931c60c7ad4SBarry Smith    Logically Collective on PC
932c60c7ad4SBarry Smith 
933c60c7ad4SBarry Smith    Input Parameter:
934c60c7ad4SBarry Smith .  pc - the preconditioner context
935c60c7ad4SBarry Smith 
936c60c7ad4SBarry Smith    Output Parameter:
937c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
938c60c7ad4SBarry Smith 
939c60c7ad4SBarry Smith 
940c60c7ad4SBarry Smith    Level: advanced
941c60c7ad4SBarry Smith 
942c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
943c60c7ad4SBarry Smith 
944c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
945c60c7ad4SBarry Smith @*/
946c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
947c60c7ad4SBarry Smith {
948c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
949c60c7ad4SBarry Smith 
950c60c7ad4SBarry Smith   PetscFunctionBegin;
951c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
952c60c7ad4SBarry Smith   *type = mg->am;
9534b9ad928SBarry Smith   PetscFunctionReturn(0);
9544b9ad928SBarry Smith }
9554b9ad928SBarry Smith 
9564b9ad928SBarry Smith /*@
9570d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
9584b9ad928SBarry Smith    complicated cycling.
9594b9ad928SBarry Smith 
960ad4df100SBarry Smith    Logically Collective on PC
9614b9ad928SBarry Smith 
9624b9ad928SBarry Smith    Input Parameters:
963c2be2410SBarry Smith +  pc - the multigrid context
964c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
9654b9ad928SBarry Smith 
9664b9ad928SBarry Smith    Options Database Key:
967c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
9684b9ad928SBarry Smith 
9694b9ad928SBarry Smith    Level: advanced
9704b9ad928SBarry Smith 
9714b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9724b9ad928SBarry Smith 
9730d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
9744b9ad928SBarry Smith @*/
9757087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
9764b9ad928SBarry Smith {
977f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
978f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
97979416396SBarry Smith   PetscInt     i,levels;
9804b9ad928SBarry Smith 
9814b9ad928SBarry Smith   PetscFunctionBegin;
9820700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
983ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
98469fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
985f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9864b9ad928SBarry Smith 
9872fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
9884b9ad928SBarry Smith   PetscFunctionReturn(0);
9894b9ad928SBarry Smith }
9904b9ad928SBarry Smith 
9918cc2d5dfSBarry Smith /*@
9928cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
9938cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9948cc2d5dfSBarry Smith 
995ad4df100SBarry Smith    Logically Collective on PC
9968cc2d5dfSBarry Smith 
9978cc2d5dfSBarry Smith    Input Parameters:
9988cc2d5dfSBarry Smith +  pc - the multigrid context
9998cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
10008cc2d5dfSBarry Smith 
10018cc2d5dfSBarry Smith    Options Database Key:
1002e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
10038cc2d5dfSBarry Smith 
10048cc2d5dfSBarry Smith    Level: advanced
10058cc2d5dfSBarry Smith 
10068cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
10078cc2d5dfSBarry Smith 
10088cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10098cc2d5dfSBarry Smith 
10108cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
10118cc2d5dfSBarry Smith @*/
10127087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
10138cc2d5dfSBarry Smith {
1014f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
10158cc2d5dfSBarry Smith 
10168cc2d5dfSBarry Smith   PetscFunctionBegin;
10170700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1018c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
10192a21e185SBarry Smith   mg->cyclesperpcapply = n;
10208cc2d5dfSBarry Smith   PetscFunctionReturn(0);
10218cc2d5dfSBarry Smith }
10228cc2d5dfSBarry Smith 
10232134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1024fb15c04eSBarry Smith {
1025fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1026fb15c04eSBarry Smith 
1027fb15c04eSBarry Smith   PetscFunctionBegin;
10282134b1e4SBarry Smith   mg->galerkin = use;
1029fb15c04eSBarry Smith   PetscFunctionReturn(0);
1030fb15c04eSBarry Smith }
1031fb15c04eSBarry Smith 
1032c2be2410SBarry Smith /*@
103397177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
103482b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1035c2be2410SBarry Smith 
1036ad4df100SBarry Smith    Logically Collective on PC
1037c2be2410SBarry Smith 
1038c2be2410SBarry Smith    Input Parameters:
1039c91913d3SJed Brown +  pc - the multigrid context
10402134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1041c2be2410SBarry Smith 
1042c2be2410SBarry Smith    Options Database Key:
10432134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1044c2be2410SBarry Smith 
1045c2be2410SBarry Smith    Level: intermediate
1046c2be2410SBarry Smith 
10471c1aac46SBarry Smith    Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
10481c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
10491c1aac46SBarry Smith 
1050c2be2410SBarry Smith .keywords: MG, set, Galerkin
1051c2be2410SBarry Smith 
10522134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
10533fc8bf9cSBarry Smith 
1054c2be2410SBarry Smith @*/
10552134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1056c2be2410SBarry Smith {
1057fb15c04eSBarry Smith   PetscErrorCode ierr;
1058c2be2410SBarry Smith 
1059c2be2410SBarry Smith   PetscFunctionBegin;
10600700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10612134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1062c2be2410SBarry Smith   PetscFunctionReturn(0);
1063c2be2410SBarry Smith }
1064c2be2410SBarry Smith 
10653fc8bf9cSBarry Smith /*@
10663fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
106782b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
10683fc8bf9cSBarry Smith 
10693fc8bf9cSBarry Smith    Not Collective
10703fc8bf9cSBarry Smith 
10713fc8bf9cSBarry Smith    Input Parameter:
10723fc8bf9cSBarry Smith .  pc - the multigrid context
10733fc8bf9cSBarry Smith 
10743fc8bf9cSBarry Smith    Output Parameter:
10752134b1e4SBarry 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
10763fc8bf9cSBarry Smith 
10773fc8bf9cSBarry Smith    Level: intermediate
10783fc8bf9cSBarry Smith 
10793fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
10803fc8bf9cSBarry Smith 
10812134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
10823fc8bf9cSBarry Smith 
10833fc8bf9cSBarry Smith @*/
10842134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
10853fc8bf9cSBarry Smith {
1086f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
10873fc8bf9cSBarry Smith 
10883fc8bf9cSBarry Smith   PetscFunctionBegin;
10890700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10902134b1e4SBarry Smith   *galerkin = mg->galerkin;
10913fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10923fc8bf9cSBarry Smith }
10933fc8bf9cSBarry Smith 
10944b9ad928SBarry Smith /*@
109597177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
109697177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
10974b9ad928SBarry Smith    pre-smoothing steps on different levels.
10984b9ad928SBarry Smith 
1099ad4df100SBarry Smith    Logically Collective on PC
11004b9ad928SBarry Smith 
11014b9ad928SBarry Smith    Input Parameters:
11024b9ad928SBarry Smith +  mg - the multigrid context
11034b9ad928SBarry Smith -  n - the number of smoothing steps
11044b9ad928SBarry Smith 
11054b9ad928SBarry Smith    Options Database Key:
11064b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
11074b9ad928SBarry Smith 
11084b9ad928SBarry Smith    Level: advanced
11098c75dd1cSBarry Smith     If the number of smoothing steps is changed in this call then the PCMGGetSmoothUp() will be called and now the
111006792cafSBarry Smith    up smoother will no longer share the same KSP object as the down smoother. Use PCMGSetNumberSmooth() to set the same
111106792cafSBarry Smith    number of smoothing steps for pre and post smoothing.
11124b9ad928SBarry Smith 
11134b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
11144b9ad928SBarry Smith 
111506792cafSBarry Smith .seealso: PCMGSetNumberSmoothUp(), PCMGSetNumberSmooth()
11164b9ad928SBarry Smith @*/
11177087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
11184b9ad928SBarry Smith {
1119f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1120f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11216849ba73SBarry Smith   PetscErrorCode ierr;
112279416396SBarry Smith   PetscInt       i,levels;
11234b9ad928SBarry Smith 
11244b9ad928SBarry Smith   PetscFunctionBegin;
11250700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1126ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1127c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1128f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11294b9ad928SBarry Smith 
1130b05257ddSBarry Smith   for (i=1; i<levels; i++) {
11318c75dd1cSBarry Smith     PetscInt nc;
11328c75dd1cSBarry Smith     ierr = KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);CHKERRQ(ierr);
11338c75dd1cSBarry Smith     if (nc == n) continue;
11348c75dd1cSBarry Smith 
11354b9ad928SBarry Smith     /* make sure smoother up and down are different */
11360298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1137f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11382fa5cd67SKarl Rupp 
113931567311SBarry Smith     mg->default_smoothd = n;
11404b9ad928SBarry Smith   }
11414b9ad928SBarry Smith   PetscFunctionReturn(0);
11424b9ad928SBarry Smith }
11434b9ad928SBarry Smith 
11444b9ad928SBarry Smith /*@
114597177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
114697177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
11474b9ad928SBarry Smith    post-smoothing steps on different levels.
11484b9ad928SBarry Smith 
1149ad4df100SBarry Smith    Logically Collective on PC
11504b9ad928SBarry Smith 
11514b9ad928SBarry Smith    Input Parameters:
11524b9ad928SBarry Smith +  mg - the multigrid context
11534b9ad928SBarry Smith -  n - the number of smoothing steps
11544b9ad928SBarry Smith 
11554b9ad928SBarry Smith    Options Database Key:
11564b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
11574b9ad928SBarry Smith 
11584b9ad928SBarry Smith    Level: advanced
11594b9ad928SBarry Smith 
11608c75dd1cSBarry Smith    Notes: this does not set a value on the coarsest grid, since we assume that
1161a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
11624b9ad928SBarry Smith 
11638c75dd1cSBarry Smith     If the number of smoothing steps is changed in this call then the PCMGGetSmoothUp() will be called and now the
116406792cafSBarry Smith    up smoother will no longer share the same KSP object as the down smoother. Use PCMGSetNumberSmooth() to set the same
116506792cafSBarry Smith    number of smoothing steps for pre and post smoothing.
116606792cafSBarry Smith 
11678c75dd1cSBarry Smith 
11684b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
11694b9ad928SBarry Smith 
117006792cafSBarry Smith .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmooth()
11714b9ad928SBarry Smith @*/
11727087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
11734b9ad928SBarry Smith {
1174f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1175f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11766849ba73SBarry Smith   PetscErrorCode ierr;
117779416396SBarry Smith   PetscInt       i,levels;
11784b9ad928SBarry Smith 
11794b9ad928SBarry Smith   PetscFunctionBegin;
11800700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1181ce94432eSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1182c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
1183f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11844b9ad928SBarry Smith 
11854b9ad928SBarry Smith   for (i=1; i<levels; i++) {
11868c75dd1cSBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
11878c75dd1cSBarry Smith       PetscInt nc;
11888c75dd1cSBarry Smith       ierr = KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);CHKERRQ(ierr);
11898c75dd1cSBarry Smith       if (nc == n) continue;
11908c75dd1cSBarry Smith     }
11918c75dd1cSBarry Smith 
11924b9ad928SBarry Smith     /* make sure smoother up and down are different */
11930298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1194f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11952fa5cd67SKarl Rupp 
119631567311SBarry Smith     mg->default_smoothu = n;
11974b9ad928SBarry Smith   }
11984b9ad928SBarry Smith   PetscFunctionReturn(0);
11994b9ad928SBarry Smith }
12004b9ad928SBarry Smith 
120106792cafSBarry Smith /*@
120206792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
120306792cafSBarry Smith    on all levels. Use PCMGSetSmoothUp() and PCMGSetSmoothDown() set different numbers of
120406792cafSBarry Smith    pre ad post-smoothing steps
120506792cafSBarry Smith 
120606792cafSBarry Smith    Logically Collective on PC
120706792cafSBarry Smith 
120806792cafSBarry Smith    Input Parameters:
120906792cafSBarry Smith +  mg - the multigrid context
121006792cafSBarry Smith -  n - the number of smoothing steps
121106792cafSBarry Smith 
121206792cafSBarry Smith    Options Database Key:
121306792cafSBarry Smith +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
121406792cafSBarry Smith .  -pc_mg_smooth_down <n> - Sets number of pre-smoothing steps (if setting different pre and post amounts)
121506792cafSBarry Smith -  -pc_mg_smooth_up <n> - Sets number of post-smoothing steps (if setting different pre and post amounts)
121606792cafSBarry Smith 
121706792cafSBarry Smith    Level: advanced
121806792cafSBarry Smith 
121906792cafSBarry Smith    Notes: this does not set a value on the coarsest grid, since we assume that
122006792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
122106792cafSBarry Smith 
122206792cafSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
122306792cafSBarry Smith 
122406792cafSBarry Smith .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmoothUp()
122506792cafSBarry Smith @*/
122606792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
122706792cafSBarry Smith {
122806792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
122906792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
123006792cafSBarry Smith   PetscErrorCode ierr;
123106792cafSBarry Smith   PetscInt       i,levels;
123206792cafSBarry Smith 
123306792cafSBarry Smith   PetscFunctionBegin;
123406792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
123506792cafSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
123606792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
123706792cafSBarry Smith   levels = mglevels[0]->levels;
123806792cafSBarry Smith 
123906792cafSBarry Smith   for (i=1; i<levels; i++) {
124006792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
124106792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
124206792cafSBarry Smith     mg->default_smoothu = n;
124306792cafSBarry Smith     mg->default_smoothd = n;
124406792cafSBarry Smith   }
124506792cafSBarry Smith   PetscFunctionReturn(0);
124606792cafSBarry Smith }
124706792cafSBarry Smith 
12484b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
12494b9ad928SBarry Smith 
12503b09bd56SBarry Smith /*MC
1251ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
12523b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
12533b09bd56SBarry Smith 
12543b09bd56SBarry Smith    Options Database Keys:
12553b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
12564afba7b0SPatrick Sanan .  -pc_mg_cycle_type <v,w> -
125779416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
12583b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
12598c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
12603b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
12612134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
12628cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
12638cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1264e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
12658cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
12668cf6037eSBarry Smith                         to the binary output file called binaryoutput
12673b09bd56SBarry Smith 
1268e941fc33SBarry Smith    Notes: If one uses a Krylov method such GMRES or CG as the smoother than one must use KSPFGMRES, KSPGCG, or KSPRICHARDSON as the outer Krylov method
12693b09bd56SBarry Smith 
12708cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
12718cf6037eSBarry Smith 
127223067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
127323067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
127423067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
127523067569SBarry Smith        residual is computed at the end of each cycle.
127623067569SBarry Smith 
12773b09bd56SBarry Smith    Level: intermediate
12783b09bd56SBarry Smith 
12798f87f92bSBarry Smith    Concepts: multigrid/multilevel
12803b09bd56SBarry Smith 
12818cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
12820d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
128397177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
128497177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12850d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12863b09bd56SBarry Smith M*/
12873b09bd56SBarry Smith 
12888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
12894b9ad928SBarry Smith {
1290f3fbd535SBarry Smith   PC_MG          *mg;
1291f3fbd535SBarry Smith   PetscErrorCode ierr;
1292f3fbd535SBarry Smith 
12934b9ad928SBarry Smith   PetscFunctionBegin;
1294b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1295f3fbd535SBarry Smith   pc->data     = (void*)mg;
1296f3fbd535SBarry Smith   mg->nlevels  = -1;
129710eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
12982134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1299f3fbd535SBarry Smith 
130037a44384SMark Adams   pc->useAmat = PETSC_TRUE;
130137a44384SMark Adams 
13024b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
13034b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1304a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
13054b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
13064b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
13074b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1308fb15c04eSBarry Smith 
1309fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
13104b9ad928SBarry Smith   PetscFunctionReturn(0);
13114b9ad928SBarry Smith }
1312