xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
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);
1360de8493bSLawrence Mitchell       ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr);
13773250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1383aeaf226SBarry Smith     }
1393aeaf226SBarry Smith 
1403aeaf226SBarry Smith     for (i=0; i<n; i++) {
1413aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1423aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1433aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1443aeaf226SBarry Smith       }
1453aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1463aeaf226SBarry Smith     }
1473aeaf226SBarry Smith   }
1483aeaf226SBarry Smith   PetscFunctionReturn(0);
1493aeaf226SBarry Smith }
1503aeaf226SBarry Smith 
1514b9ad928SBarry Smith /*@C
15297177400SBarry Smith    PCMGSetLevels - Sets the number of levels to use with MG.
1534b9ad928SBarry Smith    Must be called before any other MG routine.
1544b9ad928SBarry Smith 
155ad4df100SBarry Smith    Logically Collective on PC
1564b9ad928SBarry Smith 
1574b9ad928SBarry Smith    Input Parameters:
1584b9ad928SBarry Smith +  pc - the preconditioner context
1594b9ad928SBarry Smith .  levels - the number of levels
1604b9ad928SBarry Smith -  comms - optional communicators for each level; this is to allow solving the coarser problems
1612bf68e3eSBarry Smith            on smaller sets of processors.
1624b9ad928SBarry Smith 
1634b9ad928SBarry Smith    Level: intermediate
1644b9ad928SBarry Smith 
1654b9ad928SBarry Smith    Notes:
1664b9ad928SBarry Smith      If the number of levels is one then the multigrid uses the -mg_levels prefix
1674b9ad928SBarry Smith   for setting the level options rather than the -mg_coarse prefix.
1684b9ad928SBarry Smith 
1694b9ad928SBarry Smith .keywords: MG, set, levels, multigrid
1704b9ad928SBarry Smith 
17197177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels()
1724b9ad928SBarry Smith @*/
1737087cfbeSBarry Smith PetscErrorCode  PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
1744b9ad928SBarry Smith {
175dfbe8321SBarry Smith   PetscErrorCode ierr;
176f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
177ce94432eSBarry Smith   MPI_Comm       comm;
1783aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
17910eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
18010167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
181f3fbd535SBarry Smith   PetscInt       i;
182f3fbd535SBarry Smith   PetscMPIInt    size;
183f3fbd535SBarry Smith   const char     *prefix;
184f3fbd535SBarry Smith   PC             ipc;
1853aeaf226SBarry Smith   PetscInt       n;
1864b9ad928SBarry Smith 
1874b9ad928SBarry Smith   PetscFunctionBegin;
1880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
189548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
190548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1911a97d7b6SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1923aeaf226SBarry Smith   if (mglevels) {
19310eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
1943aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
1953aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
1963aeaf226SBarry Smith     n    = mglevels[0]->levels;
1973aeaf226SBarry Smith     for (i=0; i<n; i++) {
1983aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1993aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
2003aeaf226SBarry Smith       }
2013aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
2023aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
2033aeaf226SBarry Smith     }
2043aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
2053aeaf226SBarry Smith   }
206f3fbd535SBarry Smith 
207f3fbd535SBarry Smith   mg->nlevels = levels;
208f3fbd535SBarry Smith 
209785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
2103bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
211f3fbd535SBarry Smith 
212f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
213f3fbd535SBarry Smith 
214a9db87a2SMatthew G Knepley   mg->stageApply = 0;
215f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
216b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
2172fa5cd67SKarl Rupp 
218f3fbd535SBarry Smith     mglevels[i]->level               = i;
219f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
22010eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
221336babb1SJed Brown     mg->default_smoothu              = 2;
222336babb1SJed Brown     mg->default_smoothd              = 2;
22363e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
22463e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
22563e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
22663e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
227f3fbd535SBarry Smith 
228f3fbd535SBarry Smith     if (comms) comm = comms[i];
229f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
230422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
2310ee9a668SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
2320ee9a668SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
2330ee9a668SBarry Smith     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
2340ee9a668SBarry Smith     if (i || levels == 1) {
2350ee9a668SBarry Smith       char tprefix[128];
2360ee9a668SBarry Smith 
237336babb1SJed Brown       ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2380059c7bdSJed Brown       ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
239336babb1SJed Brown       ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
240336babb1SJed Brown       ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
241336babb1SJed Brown       ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
2420ee9a668SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
243f3fbd535SBarry Smith 
2440ee9a668SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
2450ee9a668SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
2460ee9a668SBarry Smith     } else {
247f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
248f3fbd535SBarry Smith 
2499dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
250f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
251f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
252f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
253f3fbd535SBarry Smith       if (size > 1) {
254f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
255f3fbd535SBarry Smith       } else {
256f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
257f3fbd535SBarry Smith       }
258753b7fb9SBarry Smith       ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
259f3fbd535SBarry Smith     }
2603bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
2612fa5cd67SKarl Rupp 
262f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
26331567311SBarry Smith     mg->rtol             = 0.0;
26431567311SBarry Smith     mg->abstol           = 0.0;
26531567311SBarry Smith     mg->dtol             = 0.0;
26631567311SBarry Smith     mg->ttol             = 0.0;
26731567311SBarry Smith     mg->cyclesperpcapply = 1;
268f3fbd535SBarry Smith   }
269f3fbd535SBarry Smith   mg->levels = mglevels;
27010eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2714b9ad928SBarry Smith   PetscFunctionReturn(0);
2724b9ad928SBarry Smith }
2734b9ad928SBarry Smith 
274c07bf074SBarry Smith 
275c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
276c07bf074SBarry Smith {
277c07bf074SBarry Smith   PetscErrorCode ierr;
278a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
279a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
280a06653b4SBarry Smith   PetscInt       i,n;
281c07bf074SBarry Smith 
282c07bf074SBarry Smith   PetscFunctionBegin;
283a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
284a06653b4SBarry Smith   if (mglevels) {
285a06653b4SBarry Smith     n = mglevels[0]->levels;
286a06653b4SBarry Smith     for (i=0; i<n; i++) {
287a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2886bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
289a06653b4SBarry Smith       }
2906bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
291a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
292a06653b4SBarry Smith     }
293a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
294a06653b4SBarry Smith   }
295c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
296f3fbd535SBarry Smith   PetscFunctionReturn(0);
297f3fbd535SBarry Smith }
298f3fbd535SBarry Smith 
299f3fbd535SBarry Smith 
300f3fbd535SBarry Smith 
30109573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
30209573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
30309573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
304f3fbd535SBarry Smith 
305f3fbd535SBarry Smith /*
306f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
307f3fbd535SBarry Smith              or full cycle of multigrid.
308f3fbd535SBarry Smith 
309f3fbd535SBarry Smith   Note:
310f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
311f3fbd535SBarry Smith */
312f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
313f3fbd535SBarry Smith {
314f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
315f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
316f3fbd535SBarry Smith   PetscErrorCode ierr;
317f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
318f3fbd535SBarry Smith 
319f3fbd535SBarry Smith   PetscFunctionBegin;
320a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
321e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
322298cc208SBarry Smith   for (i=0; i<levels; i++) {
323e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
32423ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
325298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
326e1d8e5deSBarry Smith     }
327e1d8e5deSBarry Smith   }
328e1d8e5deSBarry Smith 
329f3fbd535SBarry Smith   mglevels[levels-1]->b = b;
330f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
33131567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
332f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
33331567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
3340298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
335f3fbd535SBarry Smith     }
3362fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
33731567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
3382fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
33931567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
3402fa5cd67SKarl Rupp   } else {
341f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
342f3fbd535SBarry Smith   }
343a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
344f3fbd535SBarry Smith   PetscFunctionReturn(0);
345f3fbd535SBarry Smith }
346f3fbd535SBarry Smith 
347f3fbd535SBarry Smith 
3484416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
349f3fbd535SBarry Smith {
350f3fbd535SBarry Smith   PetscErrorCode   ierr;
351710315b6SLawrence Mitchell   PetscInt         levels,cycles;
3522134b1e4SBarry Smith   PetscBool        flg;
353f3fbd535SBarry Smith   PC_MG            *mg = (PC_MG*)pc->data;
3543d3eaba7SBarry Smith   PC_MG_Levels     **mglevels;
355f3fbd535SBarry Smith   PCMGType         mgtype;
356f3fbd535SBarry Smith   PCMGCycleType    mgctype;
3572134b1e4SBarry Smith   PCMGGalerkinType gtype;
358f3fbd535SBarry Smith 
359f3fbd535SBarry Smith   PetscFunctionBegin;
3601d743356SStefano Zampini   levels = PetscMax(mg->nlevels,1);
361e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
362f3fbd535SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
3631a97d7b6SStefano Zampini   if (!flg && !mg->levels && pc->dm) {
364eb3f98d2SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
365eb3f98d2SBarry Smith     levels++;
3663aeaf226SBarry Smith     mg->usedmfornumberoflevels = PETSC_TRUE;
367eb3f98d2SBarry Smith   }
3680298fd71SBarry Smith   ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
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   }
381f56b1265SBarry Smith   flg = PETSC_FALSE;
382f442ab6aSBarry Smith   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
383f442ab6aSBarry Smith   if (flg) {
384f442ab6aSBarry Smith     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
385f442ab6aSBarry Smith   }
38631567311SBarry Smith   mgtype = mg->am;
387f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
388f3fbd535SBarry Smith   if (flg) {
389f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
390f3fbd535SBarry Smith   }
39131567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
3925363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
393f3fbd535SBarry Smith     if (flg) {
394f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
395f3fbd535SBarry Smith     }
396f3fbd535SBarry Smith   }
397f3fbd535SBarry Smith   flg  = PETSC_FALSE;
3980298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
399f3fbd535SBarry Smith   if (flg) {
400f3fbd535SBarry Smith     PetscInt i;
401f3fbd535SBarry Smith     char     eventname[128];
4021a97d7b6SStefano Zampini 
403f3fbd535SBarry Smith     levels = mglevels[0]->levels;
404f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
405f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
40663e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
407f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
40863e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
409f3fbd535SBarry Smith       if (i) {
410f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
41163e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
412f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
41363e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
414f3fbd535SBarry Smith       }
415f3fbd535SBarry Smith     }
416eec35531SMatthew G Knepley 
417ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
418eec35531SMatthew G Knepley     {
419eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
420eec35531SMatthew G Knepley       PetscStageLog stageLog;
421eec35531SMatthew G Knepley       PetscInt      st;
422eec35531SMatthew G Knepley 
423eec35531SMatthew G Knepley       PetscFunctionBegin;
424eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
425eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
426eec35531SMatthew G Knepley         PetscBool same;
427eec35531SMatthew G Knepley 
428eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4292fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
430eec35531SMatthew G Knepley       }
431eec35531SMatthew G Knepley       if (!mg->stageApply) {
432eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
433eec35531SMatthew G Knepley       }
434eec35531SMatthew G Knepley     }
435ec60ca73SSatish Balay #endif
436f3fbd535SBarry Smith   }
437f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
438f3fbd535SBarry Smith   PetscFunctionReturn(0);
439f3fbd535SBarry Smith }
440f3fbd535SBarry Smith 
4416a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4426a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
4432134b1e4SBarry Smith const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0};
444f3fbd535SBarry Smith 
4459804daf3SBarry Smith #include <petscdraw.h>
446f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
447f3fbd535SBarry Smith {
448f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
449f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
450f3fbd535SBarry Smith   PetscErrorCode ierr;
451e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
452219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
453f3fbd535SBarry Smith 
454f3fbd535SBarry Smith   PetscFunctionBegin;
455251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4565b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
457219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
458f3fbd535SBarry Smith   if (iascii) {
459e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
460efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
46131567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
46231567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
463f3fbd535SBarry Smith     }
4642134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
465f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4662134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
4672134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
4682134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
4692134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
4702134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
4712134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
4724f66f45eSBarry Smith     } else {
4734f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
474f3fbd535SBarry Smith     }
4755adeb434SBarry Smith     if (mg->view){
4765adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
4775adeb434SBarry Smith     }
478f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
479f3fbd535SBarry Smith       if (!i) {
48089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
481f3fbd535SBarry Smith       } else {
48289cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
483f3fbd535SBarry Smith       }
484f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
485f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
486f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
487f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
488f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
489f3fbd535SBarry Smith       } else if (i) {
49089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
491f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
492f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
493f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
494f3fbd535SBarry Smith       }
495f3fbd535SBarry Smith     }
4965b0b0462SBarry Smith   } else if (isbinary) {
4975b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
4985b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
4995b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5005b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5015b0b0462SBarry Smith       }
5025b0b0462SBarry Smith     }
503219b1045SBarry Smith   } else if (isdraw) {
504219b1045SBarry Smith     PetscDraw draw;
505b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
506219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
507219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5080298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
509219b1045SBarry Smith     bottom = y - th;
510219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
511b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
512219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
513219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
514219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
515b4375e8dSPeter Brune       } else {
516b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
517b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
518b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
519b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
520b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
521b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
522b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
523b4375e8dSPeter Brune       }
5240298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5251cd381d6SBarry Smith       bottom -= th;
526219b1045SBarry Smith     }
5275b0b0462SBarry Smith   }
528f3fbd535SBarry Smith   PetscFunctionReturn(0);
529f3fbd535SBarry Smith }
530f3fbd535SBarry Smith 
531af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
532af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
533cab2e9ccSBarry Smith 
534f3fbd535SBarry Smith /*
535f3fbd535SBarry Smith     Calls setup for the KSP on each level
536f3fbd535SBarry Smith */
537f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
538f3fbd535SBarry Smith {
539f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
540f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
541f3fbd535SBarry Smith   PetscErrorCode ierr;
5421a97d7b6SStefano Zampini   PetscInt       i,n;
54398e04984SBarry Smith   PC             cpc;
5443db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
545f3fbd535SBarry Smith   Mat            dA,dB;
546f3fbd535SBarry Smith   Vec            tvec;
547218a07d4SBarry Smith   DM             *dms;
548649052a6SBarry Smith   PetscViewer    viewer = 0;
5492134b1e4SBarry Smith   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
550f3fbd535SBarry Smith 
551f3fbd535SBarry Smith   PetscFunctionBegin;
5521a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
5531a97d7b6SStefano Zampini   n = mglevels[0]->levels;
55401bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5553aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5563aeaf226SBarry Smith     PetscInt levels;
5573aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5583aeaf226SBarry Smith     levels++;
5593aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5600298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5613aeaf226SBarry Smith       n        = levels;
5623aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5633aeaf226SBarry Smith       mglevels =  mg->levels;
5643aeaf226SBarry Smith     }
5653aeaf226SBarry Smith   }
56698e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5673aeaf226SBarry Smith 
568f3fbd535SBarry Smith 
569f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
570f3fbd535SBarry Smith   /* so use those from global PC */
571f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
5720298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
57398e04984SBarry Smith   if (opsset) {
57498e04984SBarry Smith     Mat mmat;
57523ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
57698e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
57798e04984SBarry Smith   }
5784a5f07a5SMark F. Adams 
5794a5f07a5SMark F. Adams   if (!opsset) {
58071b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
5814a5f07a5SMark F. Adams     if(use_amat){
582f3fbd535SBarry 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);
58323ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
584f3fbd535SBarry Smith     }
5854a5f07a5SMark F. Adams     else {
5864a5f07a5SMark 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);
58723ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
5884a5f07a5SMark F. Adams     }
5894a5f07a5SMark F. Adams   }
590f3fbd535SBarry Smith 
5919c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
5929c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
5939c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
5949c7ffb39SBarry Smith       continue;
5959c7ffb39SBarry Smith     }
5969c7ffb39SBarry Smith   }
5972134b1e4SBarry Smith 
5982134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
5992134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
6002134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
6012134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
6022134b1e4SBarry Smith   }
6032134b1e4SBarry Smith 
6042134b1e4SBarry Smith 
6059c7ffb39SBarry Smith   /*
6069c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6079c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6089c7ffb39SBarry Smith   */
6092134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
6102d2b81a6SBarry Smith 	/* construct the interpolation from the DMs */
6112e499ae9SBarry Smith     Mat p;
61273250ac0SBarry Smith     Vec rscale;
613785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
614218a07d4SBarry Smith     dms[n-1] = pc->dm;
615ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
616ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
61741fce8fdSHong Zhang 	/*
61841fce8fdSHong Zhang 	   Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
61941fce8fdSHong Zhang 	   Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
62041fce8fdSHong Zhang 	   But it is safe to use -dm_mat_type.
62141fce8fdSHong Zhang 
62241fce8fdSHong Zhang 	   The mat type should not be hardcoded like this, we need to find a better way.
62341fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
62441fce8fdSHong Zhang     */
625218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
626942e3340SBarry Smith       DMKSP     kdm;
627eab52d2bSLawrence Mitchell       PetscBool dmhasrestrict, dmhasinject;
6283c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6292134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
630942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
631d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
632d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6330298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6340298fd71SBarry Smith       kdm->rhsctx          = NULL;
63524c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
636e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6372d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
63800ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
63973250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6406bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
641218a07d4SBarry Smith       }
6423ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
6433ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
6443ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
6453ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
6463ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
6473ad4599aSBarry Smith       }
648eab52d2bSLawrence Mitchell       ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr);
649eab52d2bSLawrence Mitchell       if (dmhasinject && !mglevels[i+1]->inject){
650eab52d2bSLawrence Mitchell         ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr);
651eab52d2bSLawrence Mitchell         ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr);
652eab52d2bSLawrence Mitchell         ierr = MatDestroy(&p);CHKERRQ(ierr);
653eab52d2bSLawrence Mitchell       }
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);
741eab52d2bSLawrence Mitchell       ierr   = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr);
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 
1010*95452b02SPatrick Sanan    Notes:
1011*95452b02SPatrick Sanan     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
10128cc2d5dfSBarry Smith 
10138cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10148cc2d5dfSBarry Smith 
10158cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
10168cc2d5dfSBarry Smith @*/
10177087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
10188cc2d5dfSBarry Smith {
1019f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
10208cc2d5dfSBarry Smith 
10218cc2d5dfSBarry Smith   PetscFunctionBegin;
10220700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1023c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
10242a21e185SBarry Smith   mg->cyclesperpcapply = n;
10258cc2d5dfSBarry Smith   PetscFunctionReturn(0);
10268cc2d5dfSBarry Smith }
10278cc2d5dfSBarry Smith 
10282134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1029fb15c04eSBarry Smith {
1030fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1031fb15c04eSBarry Smith 
1032fb15c04eSBarry Smith   PetscFunctionBegin;
10332134b1e4SBarry Smith   mg->galerkin = use;
1034fb15c04eSBarry Smith   PetscFunctionReturn(0);
1035fb15c04eSBarry Smith }
1036fb15c04eSBarry Smith 
1037c2be2410SBarry Smith /*@
103897177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
103982b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1040c2be2410SBarry Smith 
1041ad4df100SBarry Smith    Logically Collective on PC
1042c2be2410SBarry Smith 
1043c2be2410SBarry Smith    Input Parameters:
1044c91913d3SJed Brown +  pc - the multigrid context
10452134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1046c2be2410SBarry Smith 
1047c2be2410SBarry Smith    Options Database Key:
10482134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1049c2be2410SBarry Smith 
1050c2be2410SBarry Smith    Level: intermediate
1051c2be2410SBarry Smith 
1052*95452b02SPatrick Sanan    Notes:
1053*95452b02SPatrick Sanan     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
10541c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
10551c1aac46SBarry Smith 
1056c2be2410SBarry Smith .keywords: MG, set, Galerkin
1057c2be2410SBarry Smith 
10582134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
10593fc8bf9cSBarry Smith 
1060c2be2410SBarry Smith @*/
10612134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1062c2be2410SBarry Smith {
1063fb15c04eSBarry Smith   PetscErrorCode ierr;
1064c2be2410SBarry Smith 
1065c2be2410SBarry Smith   PetscFunctionBegin;
10660700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10672134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1068c2be2410SBarry Smith   PetscFunctionReturn(0);
1069c2be2410SBarry Smith }
1070c2be2410SBarry Smith 
10713fc8bf9cSBarry Smith /*@
10723fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
107382b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
10743fc8bf9cSBarry Smith 
10753fc8bf9cSBarry Smith    Not Collective
10763fc8bf9cSBarry Smith 
10773fc8bf9cSBarry Smith    Input Parameter:
10783fc8bf9cSBarry Smith .  pc - the multigrid context
10793fc8bf9cSBarry Smith 
10803fc8bf9cSBarry Smith    Output Parameter:
10812134b1e4SBarry 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
10823fc8bf9cSBarry Smith 
10833fc8bf9cSBarry Smith    Level: intermediate
10843fc8bf9cSBarry Smith 
10853fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
10863fc8bf9cSBarry Smith 
10872134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
10883fc8bf9cSBarry Smith 
10893fc8bf9cSBarry Smith @*/
10902134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
10913fc8bf9cSBarry Smith {
1092f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
10933fc8bf9cSBarry Smith 
10943fc8bf9cSBarry Smith   PetscFunctionBegin;
10950700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10962134b1e4SBarry Smith   *galerkin = mg->galerkin;
10973fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10983fc8bf9cSBarry Smith }
10993fc8bf9cSBarry Smith 
11004b9ad928SBarry Smith /*@
110106792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1102710315b6SLawrence Mitchell    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1103710315b6SLawrence Mitchell    pre- and post-smoothing steps.
110406792cafSBarry Smith 
110506792cafSBarry Smith    Logically Collective on PC
110606792cafSBarry Smith 
110706792cafSBarry Smith    Input Parameters:
110806792cafSBarry Smith +  mg - the multigrid context
110906792cafSBarry Smith -  n - the number of smoothing steps
111006792cafSBarry Smith 
111106792cafSBarry Smith    Options Database Key:
111206792cafSBarry Smith +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
111306792cafSBarry Smith 
111406792cafSBarry Smith    Level: advanced
111506792cafSBarry Smith 
1116*95452b02SPatrick Sanan    Notes:
1117*95452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
111806792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
111906792cafSBarry Smith 
112006792cafSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
112106792cafSBarry Smith 
1122710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp()
112306792cafSBarry Smith @*/
112406792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
112506792cafSBarry Smith {
112606792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
112706792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
112806792cafSBarry Smith   PetscErrorCode ierr;
112906792cafSBarry Smith   PetscInt       i,levels;
113006792cafSBarry Smith 
113106792cafSBarry Smith   PetscFunctionBegin;
113206792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
113306792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
11341a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
113506792cafSBarry Smith   levels = mglevels[0]->levels;
113606792cafSBarry Smith 
113706792cafSBarry Smith   for (i=1; i<levels; i++) {
113806792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
113906792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
114006792cafSBarry Smith     mg->default_smoothu = n;
114106792cafSBarry Smith     mg->default_smoothd = n;
114206792cafSBarry Smith   }
114306792cafSBarry Smith   PetscFunctionReturn(0);
114406792cafSBarry Smith }
114506792cafSBarry Smith 
1146f442ab6aSBarry Smith /*@
1147f442ab6aSBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels
1148710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1149f442ab6aSBarry Smith 
1150f442ab6aSBarry Smith    Logically Collective on PC
1151f442ab6aSBarry Smith 
1152f442ab6aSBarry Smith    Input Parameters:
1153f442ab6aSBarry Smith .  pc - the preconditioner context
1154f442ab6aSBarry Smith 
1155f442ab6aSBarry Smith    Options Database Key:
1156f442ab6aSBarry Smith .  -pc_mg_distinct_smoothup
1157f442ab6aSBarry Smith 
1158f442ab6aSBarry Smith    Level: advanced
1159f442ab6aSBarry Smith 
1160*95452b02SPatrick Sanan    Notes:
1161*95452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
1162f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1163f442ab6aSBarry Smith 
1164f442ab6aSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1165f442ab6aSBarry Smith 
1166710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth()
1167f442ab6aSBarry Smith @*/
1168f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1169f442ab6aSBarry Smith {
1170f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1171f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1172f442ab6aSBarry Smith   PetscErrorCode ierr;
1173f442ab6aSBarry Smith   PetscInt       i,levels;
1174f442ab6aSBarry Smith   KSP            subksp;
1175f442ab6aSBarry Smith 
1176f442ab6aSBarry Smith   PetscFunctionBegin;
1177f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1178f442ab6aSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1179f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1180f442ab6aSBarry Smith 
1181f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1182710315b6SLawrence Mitchell     const char *prefix = NULL;
1183f442ab6aSBarry Smith     /* make sure smoother up and down are different */
1184f442ab6aSBarry Smith     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1185710315b6SLawrence Mitchell     ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr);
1186710315b6SLawrence Mitchell     ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr);
1187f442ab6aSBarry Smith     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1188f442ab6aSBarry Smith   }
1189f442ab6aSBarry Smith   PetscFunctionReturn(0);
1190f442ab6aSBarry Smith }
1191f442ab6aSBarry Smith 
11924b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
11934b9ad928SBarry Smith 
11943b09bd56SBarry Smith /*MC
1195ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
11963b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
11973b09bd56SBarry Smith 
11983b09bd56SBarry Smith    Options Database Keys:
11993b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
12004afba7b0SPatrick Sanan .  -pc_mg_cycle_type <v,w> -
12018c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
12023b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1203710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
12042134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
12058cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
12068cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1207e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
12088cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
12098cf6037eSBarry Smith                         to the binary output file called binaryoutput
12103b09bd56SBarry Smith 
1211*95452b02SPatrick Sanan    Notes:
1212*95452b02SPatrick Sanan     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
12133b09bd56SBarry Smith 
12148cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
12158cf6037eSBarry Smith 
121623067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
121723067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
121823067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
121923067569SBarry Smith        residual is computed at the end of each cycle.
122023067569SBarry Smith 
12213b09bd56SBarry Smith    Level: intermediate
12223b09bd56SBarry Smith 
12238f87f92bSBarry Smith    Concepts: multigrid/multilevel
12243b09bd56SBarry Smith 
12258cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1226710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1227710315b6SLawrence Mitchell            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
122897177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12290d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12303b09bd56SBarry Smith M*/
12313b09bd56SBarry Smith 
12328cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
12334b9ad928SBarry Smith {
1234f3fbd535SBarry Smith   PC_MG          *mg;
1235f3fbd535SBarry Smith   PetscErrorCode ierr;
1236f3fbd535SBarry Smith 
12374b9ad928SBarry Smith   PetscFunctionBegin;
1238b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1239f3fbd535SBarry Smith   pc->data     = (void*)mg;
1240f3fbd535SBarry Smith   mg->nlevels  = -1;
124110eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
12422134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1243f3fbd535SBarry Smith 
124437a44384SMark Adams   pc->useAmat = PETSC_TRUE;
124537a44384SMark Adams 
12464b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
12474b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1248a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
12494b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
12504b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
12514b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1252fb15c04eSBarry Smith 
1253fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
12544b9ad928SBarry Smith   PetscFunctionReturn(0);
12554b9ad928SBarry Smith }
1256