xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision f442ab6a2fddf070051ab641a49a7222150491ff)
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);
189548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1901a97d7b6SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
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;
3501d743356SStefano Zampini   PetscInt         m,levels,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;
3591d743356SStefano Zampini   levels = PetscMax(mg->nlevels,1);
360e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
361f3fbd535SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
3621a97d7b6SStefano Zampini   if (!flg && !mg->levels && 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);
3683aeaf226SBarry Smith   mglevels = mg->levels;
3693aeaf226SBarry Smith 
370f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
371f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
372f3fbd535SBarry Smith   if (flg) {
373f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
3742fa5cd67SKarl Rupp   }
3752134b1e4SBarry Smith   gtype = mg->galerkin;
3762134b1e4SBarry Smith   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
3772134b1e4SBarry Smith   if (flg) {
3782134b1e4SBarry Smith     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
379f3fbd535SBarry Smith   }
38094ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",mg->default_smoothu,&m,&flg);CHKERRQ(ierr);
381f3fbd535SBarry Smith   if (flg) {
382f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr);
383f3fbd535SBarry Smith   }
38494ae4db5SBarry Smith   ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",mg->default_smoothd,&m,&flg);CHKERRQ(ierr);
385f3fbd535SBarry Smith   if (flg) {
386f3fbd535SBarry Smith     ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr);
387f3fbd535SBarry Smith   }
388*f442ab6aSBarry Smith   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
389*f442ab6aSBarry Smith   if (flg) {
390*f442ab6aSBarry Smith     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
391*f442ab6aSBarry Smith   }
39231567311SBarry Smith   mgtype = mg->am;
393f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
394f3fbd535SBarry Smith   if (flg) {
395f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
396f3fbd535SBarry Smith   }
39731567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
3985363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
399f3fbd535SBarry Smith     if (flg) {
400f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
401f3fbd535SBarry Smith     }
402f3fbd535SBarry Smith   }
403f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4040298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
405f3fbd535SBarry Smith   if (flg) {
406f3fbd535SBarry Smith     PetscInt i;
407f3fbd535SBarry Smith     char     eventname[128];
4081a97d7b6SStefano Zampini 
409f3fbd535SBarry Smith     levels = mglevels[0]->levels;
410f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
411f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
41263e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
413f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
41463e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
415f3fbd535SBarry Smith       if (i) {
416f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
41763e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
418f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
41963e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
420f3fbd535SBarry Smith       }
421f3fbd535SBarry Smith     }
422eec35531SMatthew G Knepley 
423ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
424eec35531SMatthew G Knepley     {
425eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
426eec35531SMatthew G Knepley       PetscStageLog stageLog;
427eec35531SMatthew G Knepley       PetscInt      st;
428eec35531SMatthew G Knepley 
429eec35531SMatthew G Knepley       PetscFunctionBegin;
430eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
431eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
432eec35531SMatthew G Knepley         PetscBool same;
433eec35531SMatthew G Knepley 
434eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4352fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
436eec35531SMatthew G Knepley       }
437eec35531SMatthew G Knepley       if (!mg->stageApply) {
438eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
439eec35531SMatthew G Knepley       }
440eec35531SMatthew G Knepley     }
441ec60ca73SSatish Balay #endif
442f3fbd535SBarry Smith   }
443f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
444f3fbd535SBarry Smith   PetscFunctionReturn(0);
445f3fbd535SBarry Smith }
446f3fbd535SBarry Smith 
4476a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4486a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
4492134b1e4SBarry Smith const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0};
450f3fbd535SBarry Smith 
4519804daf3SBarry Smith #include <petscdraw.h>
452f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
453f3fbd535SBarry Smith {
454f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
455f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
456f3fbd535SBarry Smith   PetscErrorCode ierr;
457e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
458219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
459f3fbd535SBarry Smith 
460f3fbd535SBarry Smith   PetscFunctionBegin;
461251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4625b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
463219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
464f3fbd535SBarry Smith   if (iascii) {
465e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
466efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
46731567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
46831567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
469f3fbd535SBarry Smith     }
4702134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
471f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4722134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
4732134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
4742134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
4752134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
4762134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
4772134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
4784f66f45eSBarry Smith     } else {
4794f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
480f3fbd535SBarry Smith     }
4815adeb434SBarry Smith     if (mg->view){
4825adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
4835adeb434SBarry Smith     }
484f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
485f3fbd535SBarry Smith       if (!i) {
48689cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
487f3fbd535SBarry Smith       } else {
48889cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
489f3fbd535SBarry Smith       }
490f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
491f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
492f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
493f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
494f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
495f3fbd535SBarry Smith       } else if (i) {
49689cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
497f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
498f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
499f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
500f3fbd535SBarry Smith       }
501f3fbd535SBarry Smith     }
5025b0b0462SBarry Smith   } else if (isbinary) {
5035b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5045b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5055b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5065b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5075b0b0462SBarry Smith       }
5085b0b0462SBarry Smith     }
509219b1045SBarry Smith   } else if (isdraw) {
510219b1045SBarry Smith     PetscDraw draw;
511b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
512219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
513219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5140298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
515219b1045SBarry Smith     bottom = y - th;
516219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
517b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
518219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
519219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
520219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
521b4375e8dSPeter Brune       } else {
522b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
523b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
524b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
525b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
526b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
527b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
528b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
529b4375e8dSPeter Brune       }
5300298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5311cd381d6SBarry Smith       bottom -= th;
532219b1045SBarry Smith     }
5335b0b0462SBarry Smith   }
534f3fbd535SBarry Smith   PetscFunctionReturn(0);
535f3fbd535SBarry Smith }
536f3fbd535SBarry Smith 
537af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
538af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
539cab2e9ccSBarry Smith 
540f3fbd535SBarry Smith /*
541f3fbd535SBarry Smith     Calls setup for the KSP on each level
542f3fbd535SBarry Smith */
543f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
544f3fbd535SBarry Smith {
545f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
546f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
547f3fbd535SBarry Smith   PetscErrorCode ierr;
5481a97d7b6SStefano Zampini   PetscInt       i,n;
54998e04984SBarry Smith   PC             cpc;
5503db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
551f3fbd535SBarry Smith   Mat            dA,dB;
552f3fbd535SBarry Smith   Vec            tvec;
553218a07d4SBarry Smith   DM             *dms;
554649052a6SBarry Smith   PetscViewer    viewer = 0;
5552134b1e4SBarry Smith   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
556f3fbd535SBarry Smith 
557f3fbd535SBarry Smith   PetscFunctionBegin;
5581a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
5591a97d7b6SStefano Zampini   n = mglevels[0]->levels;
56001bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5613aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5623aeaf226SBarry Smith     PetscInt levels;
5633aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5643aeaf226SBarry Smith     levels++;
5653aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5660298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5673aeaf226SBarry Smith       n        = levels;
5683aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5693aeaf226SBarry Smith       mglevels =  mg->levels;
5703aeaf226SBarry Smith     }
5713aeaf226SBarry Smith   }
57298e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5733aeaf226SBarry Smith 
574f3fbd535SBarry Smith 
575f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
576f3fbd535SBarry Smith   /* so use those from global PC */
577f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
5780298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
57998e04984SBarry Smith   if (opsset) {
58098e04984SBarry Smith     Mat mmat;
58123ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
58298e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
58398e04984SBarry Smith   }
5844a5f07a5SMark F. Adams 
5854a5f07a5SMark F. Adams   if (!opsset) {
58671b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
5874a5f07a5SMark F. Adams     if(use_amat){
588f3fbd535SBarry 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);
58923ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
590f3fbd535SBarry Smith     }
5914a5f07a5SMark F. Adams     else {
5924a5f07a5SMark 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);
59323ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
5944a5f07a5SMark F. Adams     }
5954a5f07a5SMark F. Adams   }
596f3fbd535SBarry Smith 
5979c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
5989c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
5999c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
6009c7ffb39SBarry Smith       continue;
6019c7ffb39SBarry Smith     }
6029c7ffb39SBarry Smith   }
6032134b1e4SBarry Smith 
6042134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
6052134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
6062134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
6072134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
6082134b1e4SBarry Smith   }
6092134b1e4SBarry Smith 
6102134b1e4SBarry Smith 
6119c7ffb39SBarry Smith   /*
6129c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6139c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6149c7ffb39SBarry Smith   */
6152134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
6162d2b81a6SBarry Smith 	/* construct the interpolation from the DMs */
6172e499ae9SBarry Smith     Mat p;
61873250ac0SBarry Smith     Vec rscale;
619785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
620218a07d4SBarry Smith     dms[n-1] = pc->dm;
621ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
622ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
62341fce8fdSHong Zhang 	/*
62441fce8fdSHong Zhang 	   Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
62541fce8fdSHong Zhang 	   Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
62641fce8fdSHong Zhang 	   But it is safe to use -dm_mat_type.
62741fce8fdSHong Zhang 
62841fce8fdSHong Zhang 	   The mat type should not be hardcoded like this, we need to find a better way.
62941fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
63041fce8fdSHong Zhang     */
631218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
632942e3340SBarry Smith       DMKSP     kdm;
6333ad4599aSBarry Smith       PetscBool dmhasrestrict;
6343c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6352134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
636942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
637d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
638d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6390298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6400298fd71SBarry Smith       kdm->rhsctx          = NULL;
64124c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
642e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6432d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
64400ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
64573250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6466bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
647218a07d4SBarry Smith       }
6483ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
6493ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
6503ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
6513ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
6523ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
6533ad4599aSBarry Smith       }
65424c3aa18SBarry Smith     }
6552d2b81a6SBarry Smith 
656ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
6572d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
658ce4cda84SJed Brown   }
659cab2e9ccSBarry Smith 
660ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
6612134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
662cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
663cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
664218a07d4SBarry Smith   }
665218a07d4SBarry Smith 
6662134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
6672134b1e4SBarry Smith     Mat       A,B;
6682134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
6692134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
6702134b1e4SBarry Smith 
6712134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
6722134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
6732134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
674f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
675b9367914SBarry 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");
676b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
677b9367914SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
678b9367914SBarry Smith       }
679b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
680b9367914SBarry Smith         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
681b9367914SBarry Smith       }
6822134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
6832134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
6842134b1e4SBarry Smith       }
6852134b1e4SBarry Smith       if (doA) {
6862df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
6872134b1e4SBarry Smith       }
6882134b1e4SBarry Smith       if (doB) {
6892df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
690b9367914SBarry Smith       }
6912134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
6922134b1e4SBarry Smith       if (!doA && dAeqdB) {
6932134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
6942134b1e4SBarry Smith         A = B;
6952134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX ) {
6962134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
6972134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
698b9367914SBarry Smith       }
6992134b1e4SBarry Smith       if (!doB && dAeqdB) {
7002134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
7012134b1e4SBarry Smith         B = A;
7022134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
70323ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
7042134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
7057d537d92SJed Brown       }
7062134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
7072134b1e4SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
7082134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
7092134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
7102134b1e4SBarry Smith       }
7112134b1e4SBarry Smith       dA = A;
712f3fbd535SBarry Smith       dB = B;
713f3fbd535SBarry Smith     }
714f3fbd535SBarry Smith   }
7152134b1e4SBarry Smith   if (needRestricts && pc->dm && pc->dm->x) {
716cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
717cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
718c88c5224SJed Brown       Mat R;
719c88c5224SJed Brown       Vec rscale;
720cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
721cab2e9ccSBarry Smith         Vec *vecs;
7222a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
723cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
724cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
725cab2e9ccSBarry Smith       }
726c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
727c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
728c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
729c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
730cab2e9ccSBarry Smith     }
731f3fbd535SBarry Smith   }
7322134b1e4SBarry Smith   if (needRestricts && pc->dm) {
733caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
734ccceb50cSJed Brown       DM  dmfine,dmcoarse;
735caa4e7f2SJed Brown       Mat Restrict,Inject;
736caa4e7f2SJed Brown       Vec rscale;
737ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
738ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
739caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
740caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
7410298fd71SBarry Smith       Inject = NULL;      /* Callback should create it if it needs Injection */
742caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
743caa4e7f2SJed Brown     }
744caa4e7f2SJed Brown   }
745f3fbd535SBarry Smith 
746f3fbd535SBarry Smith   if (!pc->setupcalled) {
747f3fbd535SBarry Smith     for (i=0; i<n; i++) {
748f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
749f3fbd535SBarry Smith     }
750f3fbd535SBarry Smith     for (i=1; i<n; i++) {
751f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
752f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
753f3fbd535SBarry Smith       }
754f3fbd535SBarry Smith     }
7553ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
756f3fbd535SBarry Smith     for (i=1; i<n; i++) {
7573ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
7583ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
759f3fbd535SBarry Smith     }
760f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
761f3fbd535SBarry Smith       if (!mglevels[i]->b) {
762f3fbd535SBarry Smith         Vec *vec;
7632a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
764f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
7656bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
766f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
767f3fbd535SBarry Smith       }
768f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
769f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
770f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
7716bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
772f3fbd535SBarry Smith       }
773f3fbd535SBarry Smith       if (!mglevels[i]->x) {
774f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
775f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
7766bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
777f3fbd535SBarry Smith       }
778f3fbd535SBarry Smith     }
779f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
780f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
781f3fbd535SBarry Smith       Vec *vec;
7822a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
783f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
7846bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
785f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
786f3fbd535SBarry Smith     }
787f3fbd535SBarry Smith   }
788f3fbd535SBarry Smith 
78998e04984SBarry Smith   if (pc->dm) {
79098e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
79198e04984SBarry Smith     for (i=0; i<n-1; i++) {
79298e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
79398e04984SBarry Smith     }
79498e04984SBarry Smith   }
795f3fbd535SBarry Smith 
796f3fbd535SBarry Smith   for (i=1; i<n; i++) {
7972a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
798f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
799f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
800f3fbd535SBarry Smith     }
80163e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
802f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
803899639b0SHong Zhang     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
804899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
805899639b0SHong Zhang     }
80663e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
807d42688cbSBarry Smith     if (!mglevels[i]->residual) {
808d42688cbSBarry Smith       Mat mat;
80913842ffbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
81054b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
811d42688cbSBarry Smith     }
812f3fbd535SBarry Smith   }
813f3fbd535SBarry Smith   for (i=1; i<n; i++) {
814f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
815f3fbd535SBarry Smith       Mat downmat,downpmat;
816f3fbd535SBarry Smith 
817f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
8180298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
819f3fbd535SBarry Smith       if (!opsset) {
82023ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
82123ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
822f3fbd535SBarry Smith       }
823f3fbd535SBarry Smith 
824f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
82563e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
826f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
827899639b0SHong Zhang       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) {
828899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
829899639b0SHong Zhang       }
83063e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
831f3fbd535SBarry Smith     }
832f3fbd535SBarry Smith   }
833f3fbd535SBarry Smith 
83463e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
835f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
836899639b0SHong Zhang   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
837899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
838899639b0SHong Zhang   }
83963e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
840f3fbd535SBarry Smith 
841f3fbd535SBarry Smith   /*
842f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
843e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
844f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
845f3fbd535SBarry Smith 
846f3fbd535SBarry Smith    Only support one or the other at the same time.
847f3fbd535SBarry Smith   */
848f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
849c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
850ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
851f3fbd535SBarry Smith   dump = PETSC_FALSE;
852f3fbd535SBarry Smith #endif
853c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
854ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
855f3fbd535SBarry Smith 
856f3fbd535SBarry Smith   if (viewer) {
857f3fbd535SBarry Smith     for (i=1; i<n; i++) {
858f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
859f3fbd535SBarry Smith     }
860f3fbd535SBarry Smith     for (i=0; i<n; i++) {
861f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
862f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
863f3fbd535SBarry Smith     }
864f3fbd535SBarry Smith   }
865f3fbd535SBarry Smith   PetscFunctionReturn(0);
866f3fbd535SBarry Smith }
867f3fbd535SBarry Smith 
868f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
869f3fbd535SBarry Smith 
8704b9ad928SBarry Smith /*@
87197177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
8724b9ad928SBarry Smith 
8734b9ad928SBarry Smith    Not Collective
8744b9ad928SBarry Smith 
8754b9ad928SBarry Smith    Input Parameter:
8764b9ad928SBarry Smith .  pc - the preconditioner context
8774b9ad928SBarry Smith 
8784b9ad928SBarry Smith    Output parameter:
8794b9ad928SBarry Smith .  levels - the number of levels
8804b9ad928SBarry Smith 
8814b9ad928SBarry Smith    Level: advanced
8824b9ad928SBarry Smith 
8834b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8844b9ad928SBarry Smith 
88597177400SBarry Smith .seealso: PCMGSetLevels()
8864b9ad928SBarry Smith @*/
8877087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8884b9ad928SBarry Smith {
889f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8904b9ad928SBarry Smith 
8914b9ad928SBarry Smith   PetscFunctionBegin;
8920700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8934482741eSBarry Smith   PetscValidIntPointer(levels,2);
894f3fbd535SBarry Smith   *levels = mg->nlevels;
8954b9ad928SBarry Smith   PetscFunctionReturn(0);
8964b9ad928SBarry Smith }
8974b9ad928SBarry Smith 
8984b9ad928SBarry Smith /*@
89997177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
9004b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
9014b9ad928SBarry Smith 
902ad4df100SBarry Smith    Logically Collective on PC
9034b9ad928SBarry Smith 
9044b9ad928SBarry Smith    Input Parameters:
9054b9ad928SBarry Smith +  pc - the preconditioner context
9069dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
9079dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
9084b9ad928SBarry Smith 
9094b9ad928SBarry Smith    Options Database Key:
9104b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
9114b9ad928SBarry Smith    additive, full, kaskade
9124b9ad928SBarry Smith 
9134b9ad928SBarry Smith    Level: advanced
9144b9ad928SBarry Smith 
9154b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
9164b9ad928SBarry Smith 
91797177400SBarry Smith .seealso: PCMGSetLevels()
9184b9ad928SBarry Smith @*/
9197087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
9204b9ad928SBarry Smith {
921f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9224b9ad928SBarry Smith 
9234b9ad928SBarry Smith   PetscFunctionBegin;
9240700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
925c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
92631567311SBarry Smith   mg->am = form;
9279dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
928c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
929c60c7ad4SBarry Smith   PetscFunctionReturn(0);
930c60c7ad4SBarry Smith }
931c60c7ad4SBarry Smith 
932c60c7ad4SBarry Smith /*@
933c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
934c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
935c60c7ad4SBarry Smith 
936c60c7ad4SBarry Smith    Logically Collective on PC
937c60c7ad4SBarry Smith 
938c60c7ad4SBarry Smith    Input Parameter:
939c60c7ad4SBarry Smith .  pc - the preconditioner context
940c60c7ad4SBarry Smith 
941c60c7ad4SBarry Smith    Output Parameter:
942c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
943c60c7ad4SBarry Smith 
944c60c7ad4SBarry Smith 
945c60c7ad4SBarry Smith    Level: advanced
946c60c7ad4SBarry Smith 
947c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
948c60c7ad4SBarry Smith 
949c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
950c60c7ad4SBarry Smith @*/
951c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
952c60c7ad4SBarry Smith {
953c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
954c60c7ad4SBarry Smith 
955c60c7ad4SBarry Smith   PetscFunctionBegin;
956c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
957c60c7ad4SBarry Smith   *type = mg->am;
9584b9ad928SBarry Smith   PetscFunctionReturn(0);
9594b9ad928SBarry Smith }
9604b9ad928SBarry Smith 
9614b9ad928SBarry Smith /*@
9620d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
9634b9ad928SBarry Smith    complicated cycling.
9644b9ad928SBarry Smith 
965ad4df100SBarry Smith    Logically Collective on PC
9664b9ad928SBarry Smith 
9674b9ad928SBarry Smith    Input Parameters:
968c2be2410SBarry Smith +  pc - the multigrid context
969c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
9704b9ad928SBarry Smith 
9714b9ad928SBarry Smith    Options Database Key:
972c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
9734b9ad928SBarry Smith 
9744b9ad928SBarry Smith    Level: advanced
9754b9ad928SBarry Smith 
9764b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9774b9ad928SBarry Smith 
9780d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
9794b9ad928SBarry Smith @*/
9807087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
9814b9ad928SBarry Smith {
982f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
983f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
98479416396SBarry Smith   PetscInt     i,levels;
9854b9ad928SBarry Smith 
9864b9ad928SBarry Smith   PetscFunctionBegin;
9870700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
98869fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
9891a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
990f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9912fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
9924b9ad928SBarry Smith   PetscFunctionReturn(0);
9934b9ad928SBarry Smith }
9944b9ad928SBarry Smith 
9958cc2d5dfSBarry Smith /*@
9968cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
9978cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9988cc2d5dfSBarry Smith 
999ad4df100SBarry Smith    Logically Collective on PC
10008cc2d5dfSBarry Smith 
10018cc2d5dfSBarry Smith    Input Parameters:
10028cc2d5dfSBarry Smith +  pc - the multigrid context
10038cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
10048cc2d5dfSBarry Smith 
10058cc2d5dfSBarry Smith    Options Database Key:
1006e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
10078cc2d5dfSBarry Smith 
10088cc2d5dfSBarry Smith    Level: advanced
10098cc2d5dfSBarry Smith 
10108cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
10118cc2d5dfSBarry Smith 
10128cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10138cc2d5dfSBarry Smith 
10148cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
10158cc2d5dfSBarry Smith @*/
10167087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
10178cc2d5dfSBarry Smith {
1018f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
10198cc2d5dfSBarry Smith 
10208cc2d5dfSBarry Smith   PetscFunctionBegin;
10210700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1022c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
10232a21e185SBarry Smith   mg->cyclesperpcapply = n;
10248cc2d5dfSBarry Smith   PetscFunctionReturn(0);
10258cc2d5dfSBarry Smith }
10268cc2d5dfSBarry Smith 
10272134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1028fb15c04eSBarry Smith {
1029fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1030fb15c04eSBarry Smith 
1031fb15c04eSBarry Smith   PetscFunctionBegin;
10322134b1e4SBarry Smith   mg->galerkin = use;
1033fb15c04eSBarry Smith   PetscFunctionReturn(0);
1034fb15c04eSBarry Smith }
1035fb15c04eSBarry Smith 
1036c2be2410SBarry Smith /*@
103797177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
103882b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1039c2be2410SBarry Smith 
1040ad4df100SBarry Smith    Logically Collective on PC
1041c2be2410SBarry Smith 
1042c2be2410SBarry Smith    Input Parameters:
1043c91913d3SJed Brown +  pc - the multigrid context
10442134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1045c2be2410SBarry Smith 
1046c2be2410SBarry Smith    Options Database Key:
10472134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1048c2be2410SBarry Smith 
1049c2be2410SBarry Smith    Level: intermediate
1050c2be2410SBarry Smith 
10511c1aac46SBarry Smith    Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
10521c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
10531c1aac46SBarry Smith 
1054c2be2410SBarry Smith .keywords: MG, set, Galerkin
1055c2be2410SBarry Smith 
10562134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
10573fc8bf9cSBarry Smith 
1058c2be2410SBarry Smith @*/
10592134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1060c2be2410SBarry Smith {
1061fb15c04eSBarry Smith   PetscErrorCode ierr;
1062c2be2410SBarry Smith 
1063c2be2410SBarry Smith   PetscFunctionBegin;
10640700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10652134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1066c2be2410SBarry Smith   PetscFunctionReturn(0);
1067c2be2410SBarry Smith }
1068c2be2410SBarry Smith 
10693fc8bf9cSBarry Smith /*@
10703fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
107182b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
10723fc8bf9cSBarry Smith 
10733fc8bf9cSBarry Smith    Not Collective
10743fc8bf9cSBarry Smith 
10753fc8bf9cSBarry Smith    Input Parameter:
10763fc8bf9cSBarry Smith .  pc - the multigrid context
10773fc8bf9cSBarry Smith 
10783fc8bf9cSBarry Smith    Output Parameter:
10792134b1e4SBarry 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
10803fc8bf9cSBarry Smith 
10813fc8bf9cSBarry Smith    Level: intermediate
10823fc8bf9cSBarry Smith 
10833fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
10843fc8bf9cSBarry Smith 
10852134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
10863fc8bf9cSBarry Smith 
10873fc8bf9cSBarry Smith @*/
10882134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
10893fc8bf9cSBarry Smith {
1090f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
10913fc8bf9cSBarry Smith 
10923fc8bf9cSBarry Smith   PetscFunctionBegin;
10930700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10942134b1e4SBarry Smith   *galerkin = mg->galerkin;
10953fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10963fc8bf9cSBarry Smith }
10973fc8bf9cSBarry Smith 
10984b9ad928SBarry Smith /*@
109997177400SBarry Smith    PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
110097177400SBarry Smith    use on all levels. Use PCMGGetSmootherDown() to set different
11014b9ad928SBarry Smith    pre-smoothing steps on different levels.
11024b9ad928SBarry Smith 
1103ad4df100SBarry Smith    Logically Collective on PC
11044b9ad928SBarry Smith 
11054b9ad928SBarry Smith    Input Parameters:
11064b9ad928SBarry Smith +  mg - the multigrid context
11074b9ad928SBarry Smith -  n - the number of smoothing steps
11084b9ad928SBarry Smith 
11094b9ad928SBarry Smith    Options Database Key:
11104b9ad928SBarry Smith .  -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
11114b9ad928SBarry Smith 
11124b9ad928SBarry Smith    Level: advanced
11138c75dd1cSBarry Smith     If the number of smoothing steps is changed in this call then the PCMGGetSmoothUp() will be called and now the
111406792cafSBarry Smith    up smoother will no longer share the same KSP object as the down smoother. Use PCMGSetNumberSmooth() to set the same
111506792cafSBarry Smith    number of smoothing steps for pre and post smoothing.
11164b9ad928SBarry Smith 
11174b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
11184b9ad928SBarry Smith 
111906792cafSBarry Smith .seealso: PCMGSetNumberSmoothUp(), PCMGSetNumberSmooth()
11204b9ad928SBarry Smith @*/
11217087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothDown(PC pc,PetscInt n)
11224b9ad928SBarry Smith {
1123f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1124f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11256849ba73SBarry Smith   PetscErrorCode ierr;
112679416396SBarry Smith   PetscInt       i,levels;
11274b9ad928SBarry Smith 
11284b9ad928SBarry Smith   PetscFunctionBegin;
11290700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1130c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
11311a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1132f3fbd535SBarry Smith   levels = mglevels[0]->levels;
1133b05257ddSBarry Smith   for (i=1; i<levels; i++) {
11348c75dd1cSBarry Smith     PetscInt nc;
11358c75dd1cSBarry Smith     ierr = KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);CHKERRQ(ierr);
11368c75dd1cSBarry Smith     if (nc == n) continue;
11378c75dd1cSBarry Smith 
11384b9ad928SBarry Smith     /* make sure smoother up and down are different */
11390298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1140f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11412fa5cd67SKarl Rupp 
114231567311SBarry Smith     mg->default_smoothd = n;
11434b9ad928SBarry Smith   }
11444b9ad928SBarry Smith   PetscFunctionReturn(0);
11454b9ad928SBarry Smith }
11464b9ad928SBarry Smith 
11474b9ad928SBarry Smith /*@
114897177400SBarry Smith    PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
114997177400SBarry Smith    on all levels. Use PCMGGetSmootherUp() to set different numbers of
11504b9ad928SBarry Smith    post-smoothing steps on different levels.
11514b9ad928SBarry Smith 
1152ad4df100SBarry Smith    Logically Collective on PC
11534b9ad928SBarry Smith 
11544b9ad928SBarry Smith    Input Parameters:
11554b9ad928SBarry Smith +  mg - the multigrid context
11564b9ad928SBarry Smith -  n - the number of smoothing steps
11574b9ad928SBarry Smith 
11584b9ad928SBarry Smith    Options Database Key:
11594b9ad928SBarry Smith .  -pc_mg_smoothup <n> - Sets number of post-smoothing steps
11604b9ad928SBarry Smith 
11614b9ad928SBarry Smith    Level: advanced
11624b9ad928SBarry Smith 
11638c75dd1cSBarry Smith    Notes: this does not set a value on the coarsest grid, since we assume that
1164a8c7a070SBarry Smith     there is no separate smooth up on the coarsest grid.
11654b9ad928SBarry Smith 
11668c75dd1cSBarry Smith     If the number of smoothing steps is changed in this call then the PCMGGetSmoothUp() will be called and now the
116706792cafSBarry Smith    up smoother will no longer share the same KSP object as the down smoother. Use PCMGSetNumberSmooth() to set the same
116806792cafSBarry Smith    number of smoothing steps for pre and post smoothing.
116906792cafSBarry Smith 
11708c75dd1cSBarry Smith 
11714b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
11724b9ad928SBarry Smith 
117306792cafSBarry Smith .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmooth()
11744b9ad928SBarry Smith @*/
11757087cfbeSBarry Smith PetscErrorCode  PCMGSetNumberSmoothUp(PC pc,PetscInt n)
11764b9ad928SBarry Smith {
1177f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1178f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
11796849ba73SBarry Smith   PetscErrorCode ierr;
118079416396SBarry Smith   PetscInt       i,levels;
11814b9ad928SBarry Smith 
11824b9ad928SBarry Smith   PetscFunctionBegin;
11830700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1184c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
11851a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1186f3fbd535SBarry Smith   levels = mglevels[0]->levels;
11874b9ad928SBarry Smith 
11884b9ad928SBarry Smith   for (i=1; i<levels; i++) {
11898c75dd1cSBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
11908c75dd1cSBarry Smith       PetscInt nc;
11918c75dd1cSBarry Smith       ierr = KSPGetTolerances(mglevels[i]->smoothd,NULL,NULL,NULL,&nc);CHKERRQ(ierr);
11928c75dd1cSBarry Smith       if (nc == n) continue;
11938c75dd1cSBarry Smith     }
11948c75dd1cSBarry Smith 
11954b9ad928SBarry Smith     /* make sure smoother up and down are different */
11960298fd71SBarry Smith     ierr = PCMGGetSmootherUp(pc,i,NULL);CHKERRQ(ierr);
1197f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
11982fa5cd67SKarl Rupp 
119931567311SBarry Smith     mg->default_smoothu = n;
12004b9ad928SBarry Smith   }
12014b9ad928SBarry Smith   PetscFunctionReturn(0);
12024b9ad928SBarry Smith }
12034b9ad928SBarry Smith 
120406792cafSBarry Smith /*@
120506792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
120606792cafSBarry Smith    on all levels. Use PCMGSetSmoothUp() and PCMGSetSmoothDown() set different numbers of
120706792cafSBarry Smith    pre ad post-smoothing steps
120806792cafSBarry Smith 
120906792cafSBarry Smith    Logically Collective on PC
121006792cafSBarry Smith 
121106792cafSBarry Smith    Input Parameters:
121206792cafSBarry Smith +  mg - the multigrid context
121306792cafSBarry Smith -  n - the number of smoothing steps
121406792cafSBarry Smith 
121506792cafSBarry Smith    Options Database Key:
121606792cafSBarry Smith +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
121706792cafSBarry Smith .  -pc_mg_smooth_down <n> - Sets number of pre-smoothing steps (if setting different pre and post amounts)
121806792cafSBarry Smith -  -pc_mg_smooth_up <n> - Sets number of post-smoothing steps (if setting different pre and post amounts)
121906792cafSBarry Smith 
122006792cafSBarry Smith    Level: advanced
122106792cafSBarry Smith 
122206792cafSBarry Smith    Notes: this does not set a value on the coarsest grid, since we assume that
122306792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
122406792cafSBarry Smith 
122506792cafSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
122606792cafSBarry Smith 
122706792cafSBarry Smith .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmoothUp()
122806792cafSBarry Smith @*/
122906792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
123006792cafSBarry Smith {
123106792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
123206792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
123306792cafSBarry Smith   PetscErrorCode ierr;
123406792cafSBarry Smith   PetscInt       i,levels;
123506792cafSBarry Smith 
123606792cafSBarry Smith   PetscFunctionBegin;
123706792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
123806792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
12391a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
124006792cafSBarry Smith   levels = mglevels[0]->levels;
124106792cafSBarry Smith 
124206792cafSBarry Smith   for (i=1; i<levels; i++) {
124306792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
124406792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
124506792cafSBarry Smith     mg->default_smoothu = n;
124606792cafSBarry Smith     mg->default_smoothd = n;
124706792cafSBarry Smith   }
124806792cafSBarry Smith   PetscFunctionReturn(0);
124906792cafSBarry Smith }
125006792cafSBarry Smith 
1251*f442ab6aSBarry Smith /*@
1252*f442ab6aSBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels
1253*f442ab6aSBarry Smith        and adds the suffix _post to the options name
1254*f442ab6aSBarry Smith 
1255*f442ab6aSBarry Smith    Logically Collective on PC
1256*f442ab6aSBarry Smith 
1257*f442ab6aSBarry Smith    Input Parameters:
1258*f442ab6aSBarry Smith .  pc - the preconditioner context
1259*f442ab6aSBarry Smith 
1260*f442ab6aSBarry Smith    Options Database Key:
1261*f442ab6aSBarry Smith .  -pc_mg_distinct_smoothup
1262*f442ab6aSBarry Smith 
1263*f442ab6aSBarry Smith    Level: advanced
1264*f442ab6aSBarry Smith 
1265*f442ab6aSBarry Smith    Notes: this does not set a value on the coarsest grid, since we assume that
1266*f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1267*f442ab6aSBarry Smith 
1268*f442ab6aSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1269*f442ab6aSBarry Smith 
1270*f442ab6aSBarry Smith .seealso: PCMGSetNumberSmoothDown(), PCMGSetNumberSmooth()
1271*f442ab6aSBarry Smith @*/
1272*f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1273*f442ab6aSBarry Smith {
1274*f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1275*f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1276*f442ab6aSBarry Smith   PetscErrorCode ierr;
1277*f442ab6aSBarry Smith   PetscInt       i,levels;
1278*f442ab6aSBarry Smith   KSP            subksp;
1279*f442ab6aSBarry Smith 
1280*f442ab6aSBarry Smith   PetscFunctionBegin;
1281*f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1282*f442ab6aSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1283*f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1284*f442ab6aSBarry Smith 
1285*f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1286*f442ab6aSBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) continue;
1287*f442ab6aSBarry Smith 
1288*f442ab6aSBarry Smith     /* make sure smoother up and down are different */
1289*f442ab6aSBarry Smith     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1290*f442ab6aSBarry Smith     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1291*f442ab6aSBarry Smith   }
1292*f442ab6aSBarry Smith   PetscFunctionReturn(0);
1293*f442ab6aSBarry Smith }
1294*f442ab6aSBarry Smith 
12954b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
12964b9ad928SBarry Smith 
12973b09bd56SBarry Smith /*MC
1298ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
12993b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
13003b09bd56SBarry Smith 
13013b09bd56SBarry Smith    Options Database Keys:
13023b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
13034afba7b0SPatrick Sanan .  -pc_mg_cycle_type <v,w> -
130479416396SBarry Smith .  -pc_mg_smoothup <n> - number of smoothing steps after interpolation
13053b09bd56SBarry Smith .  -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
13068c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
13073b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
13082134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
13098cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
13108cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1311e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
13128cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
13138cf6037eSBarry Smith                         to the binary output file called binaryoutput
13143b09bd56SBarry Smith 
1315e941fc33SBarry 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
13163b09bd56SBarry Smith 
13178cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
13188cf6037eSBarry Smith 
131923067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
132023067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
132123067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
132223067569SBarry Smith        residual is computed at the end of each cycle.
132323067569SBarry Smith 
13243b09bd56SBarry Smith    Level: intermediate
13253b09bd56SBarry Smith 
13268f87f92bSBarry Smith    Concepts: multigrid/multilevel
13273b09bd56SBarry Smith 
13288cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
13290d353602SBarry Smith            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
133097177400SBarry Smith            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
133197177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
13320d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
13333b09bd56SBarry Smith M*/
13343b09bd56SBarry Smith 
13358cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
13364b9ad928SBarry Smith {
1337f3fbd535SBarry Smith   PC_MG          *mg;
1338f3fbd535SBarry Smith   PetscErrorCode ierr;
1339f3fbd535SBarry Smith 
13404b9ad928SBarry Smith   PetscFunctionBegin;
1341b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1342f3fbd535SBarry Smith   pc->data     = (void*)mg;
1343f3fbd535SBarry Smith   mg->nlevels  = -1;
134410eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
13452134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1346f3fbd535SBarry Smith 
134737a44384SMark Adams   pc->useAmat = PETSC_TRUE;
134837a44384SMark Adams 
13494b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
13504b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1351a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
13524b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
13534b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
13544b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1355fb15c04eSBarry Smith 
1356fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
13574b9ad928SBarry Smith   PetscFunctionReturn(0);
13584b9ad928SBarry Smith }
1359