xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 710315b6501bf79f707912fe301a8dc5ebd0cb6d)
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;
350*710315b6SLawrence Mitchell   PetscInt         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   }
380f442ab6aSBarry Smith   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
381f442ab6aSBarry Smith   if (flg) {
382f442ab6aSBarry Smith     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
383f442ab6aSBarry Smith   }
38431567311SBarry Smith   mgtype = mg->am;
385f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
386f3fbd535SBarry Smith   if (flg) {
387f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
388f3fbd535SBarry Smith   }
38931567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
3905363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
391f3fbd535SBarry Smith     if (flg) {
392f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
393f3fbd535SBarry Smith     }
394f3fbd535SBarry Smith   }
395f3fbd535SBarry Smith   flg  = PETSC_FALSE;
3960298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
397f3fbd535SBarry Smith   if (flg) {
398f3fbd535SBarry Smith     PetscInt i;
399f3fbd535SBarry Smith     char     eventname[128];
4001a97d7b6SStefano Zampini 
401f3fbd535SBarry Smith     levels = mglevels[0]->levels;
402f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
403f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
40463e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
405f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
40663e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
407f3fbd535SBarry Smith       if (i) {
408f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
40963e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
410f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
41163e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
412f3fbd535SBarry Smith       }
413f3fbd535SBarry Smith     }
414eec35531SMatthew G Knepley 
415ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
416eec35531SMatthew G Knepley     {
417eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
418eec35531SMatthew G Knepley       PetscStageLog stageLog;
419eec35531SMatthew G Knepley       PetscInt      st;
420eec35531SMatthew G Knepley 
421eec35531SMatthew G Knepley       PetscFunctionBegin;
422eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
423eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
424eec35531SMatthew G Knepley         PetscBool same;
425eec35531SMatthew G Knepley 
426eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4272fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
428eec35531SMatthew G Knepley       }
429eec35531SMatthew G Knepley       if (!mg->stageApply) {
430eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
431eec35531SMatthew G Knepley       }
432eec35531SMatthew G Knepley     }
433ec60ca73SSatish Balay #endif
434f3fbd535SBarry Smith   }
435f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
436f3fbd535SBarry Smith   PetscFunctionReturn(0);
437f3fbd535SBarry Smith }
438f3fbd535SBarry Smith 
4396a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4406a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
4412134b1e4SBarry Smith const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0};
442f3fbd535SBarry Smith 
4439804daf3SBarry Smith #include <petscdraw.h>
444f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
445f3fbd535SBarry Smith {
446f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
447f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
448f3fbd535SBarry Smith   PetscErrorCode ierr;
449e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
450219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
451f3fbd535SBarry Smith 
452f3fbd535SBarry Smith   PetscFunctionBegin;
453251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
4545b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
455219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
456f3fbd535SBarry Smith   if (iascii) {
457e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
458efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
45931567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
46031567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
461f3fbd535SBarry Smith     }
4622134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
463f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
4642134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
4652134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
4662134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
4672134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
4682134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
4692134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
4704f66f45eSBarry Smith     } else {
4714f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
472f3fbd535SBarry Smith     }
4735adeb434SBarry Smith     if (mg->view){
4745adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
4755adeb434SBarry Smith     }
476f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
477f3fbd535SBarry Smith       if (!i) {
47889cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
479f3fbd535SBarry Smith       } else {
48089cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
481f3fbd535SBarry Smith       }
482f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
483f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
484f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
485f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
486f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
487f3fbd535SBarry Smith       } else if (i) {
48889cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
489f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
490f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
491f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
492f3fbd535SBarry Smith       }
493f3fbd535SBarry Smith     }
4945b0b0462SBarry Smith   } else if (isbinary) {
4955b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
4965b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
4975b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
4985b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
4995b0b0462SBarry Smith       }
5005b0b0462SBarry Smith     }
501219b1045SBarry Smith   } else if (isdraw) {
502219b1045SBarry Smith     PetscDraw draw;
503b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
504219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
505219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5060298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
507219b1045SBarry Smith     bottom = y - th;
508219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
509b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
510219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
511219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
512219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
513b4375e8dSPeter Brune       } else {
514b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
515b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
516b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
517b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
518b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
519b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
520b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
521b4375e8dSPeter Brune       }
5220298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5231cd381d6SBarry Smith       bottom -= th;
524219b1045SBarry Smith     }
5255b0b0462SBarry Smith   }
526f3fbd535SBarry Smith   PetscFunctionReturn(0);
527f3fbd535SBarry Smith }
528f3fbd535SBarry Smith 
529af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
530af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
531cab2e9ccSBarry Smith 
532f3fbd535SBarry Smith /*
533f3fbd535SBarry Smith     Calls setup for the KSP on each level
534f3fbd535SBarry Smith */
535f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
536f3fbd535SBarry Smith {
537f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
538f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
539f3fbd535SBarry Smith   PetscErrorCode ierr;
5401a97d7b6SStefano Zampini   PetscInt       i,n;
54198e04984SBarry Smith   PC             cpc;
5423db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
543f3fbd535SBarry Smith   Mat            dA,dB;
544f3fbd535SBarry Smith   Vec            tvec;
545218a07d4SBarry Smith   DM             *dms;
546649052a6SBarry Smith   PetscViewer    viewer = 0;
5472134b1e4SBarry Smith   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
548f3fbd535SBarry Smith 
549f3fbd535SBarry Smith   PetscFunctionBegin;
5501a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
5511a97d7b6SStefano Zampini   n = mglevels[0]->levels;
55201bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
5533aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
5543aeaf226SBarry Smith     PetscInt levels;
5553aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
5563aeaf226SBarry Smith     levels++;
5573aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
5580298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
5593aeaf226SBarry Smith       n        = levels;
5603aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
5613aeaf226SBarry Smith       mglevels =  mg->levels;
5623aeaf226SBarry Smith     }
5633aeaf226SBarry Smith   }
56498e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
5653aeaf226SBarry Smith 
566f3fbd535SBarry Smith 
567f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
568f3fbd535SBarry Smith   /* so use those from global PC */
569f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
5700298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
57198e04984SBarry Smith   if (opsset) {
57298e04984SBarry Smith     Mat mmat;
57323ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
57498e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
57598e04984SBarry Smith   }
5764a5f07a5SMark F. Adams 
5774a5f07a5SMark F. Adams   if (!opsset) {
57871b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
5794a5f07a5SMark F. Adams     if(use_amat){
580f3fbd535SBarry 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);
58123ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
582f3fbd535SBarry Smith     }
5834a5f07a5SMark F. Adams     else {
5844a5f07a5SMark 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);
58523ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
5864a5f07a5SMark F. Adams     }
5874a5f07a5SMark F. Adams   }
588f3fbd535SBarry Smith 
5899c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
5909c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
5919c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
5929c7ffb39SBarry Smith       continue;
5939c7ffb39SBarry Smith     }
5949c7ffb39SBarry Smith   }
5952134b1e4SBarry Smith 
5962134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
5972134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
5982134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
5992134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
6002134b1e4SBarry Smith   }
6012134b1e4SBarry Smith 
6022134b1e4SBarry Smith 
6039c7ffb39SBarry Smith   /*
6049c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6059c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6069c7ffb39SBarry Smith   */
6072134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
6082d2b81a6SBarry Smith 	/* construct the interpolation from the DMs */
6092e499ae9SBarry Smith     Mat p;
61073250ac0SBarry Smith     Vec rscale;
611785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
612218a07d4SBarry Smith     dms[n-1] = pc->dm;
613ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
614ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
61541fce8fdSHong Zhang 	/*
61641fce8fdSHong Zhang 	   Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
61741fce8fdSHong Zhang 	   Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
61841fce8fdSHong Zhang 	   But it is safe to use -dm_mat_type.
61941fce8fdSHong Zhang 
62041fce8fdSHong Zhang 	   The mat type should not be hardcoded like this, we need to find a better way.
62141fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
62241fce8fdSHong Zhang     */
623218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
624942e3340SBarry Smith       DMKSP     kdm;
6253ad4599aSBarry Smith       PetscBool dmhasrestrict;
6263c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6272134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
628942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
629d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
630d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6310298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6320298fd71SBarry Smith       kdm->rhsctx          = NULL;
63324c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
634e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6352d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
63600ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
63773250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6386bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
639218a07d4SBarry Smith       }
6403ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
6413ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
6423ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
6433ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
6443ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
6453ad4599aSBarry Smith       }
64624c3aa18SBarry Smith     }
6472d2b81a6SBarry Smith 
648ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
6492d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
650ce4cda84SJed Brown   }
651cab2e9ccSBarry Smith 
652ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
6532134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
654cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
655cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
656218a07d4SBarry Smith   }
657218a07d4SBarry Smith 
6582134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
6592134b1e4SBarry Smith     Mat       A,B;
6602134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
6612134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
6622134b1e4SBarry Smith 
6632134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
6642134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
6652134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
666f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
667b9367914SBarry 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");
668b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
669b9367914SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
670b9367914SBarry Smith       }
671b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
672b9367914SBarry Smith         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
673b9367914SBarry Smith       }
6742134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
6752134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
6762134b1e4SBarry Smith       }
6772134b1e4SBarry Smith       if (doA) {
6782df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
6792134b1e4SBarry Smith       }
6802134b1e4SBarry Smith       if (doB) {
6812df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
682b9367914SBarry Smith       }
6832134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
6842134b1e4SBarry Smith       if (!doA && dAeqdB) {
6852134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
6862134b1e4SBarry Smith         A = B;
6872134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX ) {
6882134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
6892134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
690b9367914SBarry Smith       }
6912134b1e4SBarry Smith       if (!doB && dAeqdB) {
6922134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
6932134b1e4SBarry Smith         B = A;
6942134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
69523ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
6962134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
6977d537d92SJed Brown       }
6982134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
6992134b1e4SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
7002134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
7012134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
7022134b1e4SBarry Smith       }
7032134b1e4SBarry Smith       dA = A;
704f3fbd535SBarry Smith       dB = B;
705f3fbd535SBarry Smith     }
706f3fbd535SBarry Smith   }
7072134b1e4SBarry Smith   if (needRestricts && pc->dm && pc->dm->x) {
708cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
709cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
710c88c5224SJed Brown       Mat R;
711c88c5224SJed Brown       Vec rscale;
712cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
713cab2e9ccSBarry Smith         Vec *vecs;
7142a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
715cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
716cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
717cab2e9ccSBarry Smith       }
718c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
719c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
720c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
721c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
722cab2e9ccSBarry Smith     }
723f3fbd535SBarry Smith   }
7242134b1e4SBarry Smith   if (needRestricts && pc->dm) {
725caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
726ccceb50cSJed Brown       DM  dmfine,dmcoarse;
727caa4e7f2SJed Brown       Mat Restrict,Inject;
728caa4e7f2SJed Brown       Vec rscale;
729ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
730ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
731caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
732caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
7330298fd71SBarry Smith       Inject = NULL;      /* Callback should create it if it needs Injection */
734caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
735caa4e7f2SJed Brown     }
736caa4e7f2SJed Brown   }
737f3fbd535SBarry Smith 
738f3fbd535SBarry Smith   if (!pc->setupcalled) {
739f3fbd535SBarry Smith     for (i=0; i<n; i++) {
740f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
741f3fbd535SBarry Smith     }
742f3fbd535SBarry Smith     for (i=1; i<n; i++) {
743f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
744f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
745f3fbd535SBarry Smith       }
746f3fbd535SBarry Smith     }
7473ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
748f3fbd535SBarry Smith     for (i=1; i<n; i++) {
7493ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
7503ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
751f3fbd535SBarry Smith     }
752f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
753f3fbd535SBarry Smith       if (!mglevels[i]->b) {
754f3fbd535SBarry Smith         Vec *vec;
7552a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
756f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
7576bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
758f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
759f3fbd535SBarry Smith       }
760f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
761f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
762f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
7636bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
764f3fbd535SBarry Smith       }
765f3fbd535SBarry Smith       if (!mglevels[i]->x) {
766f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
767f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
7686bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
769f3fbd535SBarry Smith       }
770f3fbd535SBarry Smith     }
771f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
772f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
773f3fbd535SBarry Smith       Vec *vec;
7742a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
775f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
7766bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
777f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
778f3fbd535SBarry Smith     }
779f3fbd535SBarry Smith   }
780f3fbd535SBarry Smith 
78198e04984SBarry Smith   if (pc->dm) {
78298e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
78398e04984SBarry Smith     for (i=0; i<n-1; i++) {
78498e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
78598e04984SBarry Smith     }
78698e04984SBarry Smith   }
787f3fbd535SBarry Smith 
788f3fbd535SBarry Smith   for (i=1; i<n; i++) {
7892a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
790f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
791f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
792f3fbd535SBarry Smith     }
79363e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
794f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
795899639b0SHong Zhang     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
796899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
797899639b0SHong Zhang     }
79863e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
799d42688cbSBarry Smith     if (!mglevels[i]->residual) {
800d42688cbSBarry Smith       Mat mat;
80113842ffbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
80254b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
803d42688cbSBarry Smith     }
804f3fbd535SBarry Smith   }
805f3fbd535SBarry Smith   for (i=1; i<n; i++) {
806f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
807f3fbd535SBarry Smith       Mat downmat,downpmat;
808f3fbd535SBarry Smith 
809f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
8100298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
811f3fbd535SBarry Smith       if (!opsset) {
81223ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
81323ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
814f3fbd535SBarry Smith       }
815f3fbd535SBarry Smith 
816f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
81763e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
818f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
819899639b0SHong Zhang       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) {
820899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
821899639b0SHong Zhang       }
82263e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
823f3fbd535SBarry Smith     }
824f3fbd535SBarry Smith   }
825f3fbd535SBarry Smith 
82663e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
827f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
828899639b0SHong Zhang   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
829899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
830899639b0SHong Zhang   }
83163e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
832f3fbd535SBarry Smith 
833f3fbd535SBarry Smith   /*
834f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
835e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
836f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
837f3fbd535SBarry Smith 
838f3fbd535SBarry Smith    Only support one or the other at the same time.
839f3fbd535SBarry Smith   */
840f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
841c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
842ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
843f3fbd535SBarry Smith   dump = PETSC_FALSE;
844f3fbd535SBarry Smith #endif
845c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
846ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
847f3fbd535SBarry Smith 
848f3fbd535SBarry Smith   if (viewer) {
849f3fbd535SBarry Smith     for (i=1; i<n; i++) {
850f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
851f3fbd535SBarry Smith     }
852f3fbd535SBarry Smith     for (i=0; i<n; i++) {
853f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
854f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
855f3fbd535SBarry Smith     }
856f3fbd535SBarry Smith   }
857f3fbd535SBarry Smith   PetscFunctionReturn(0);
858f3fbd535SBarry Smith }
859f3fbd535SBarry Smith 
860f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
861f3fbd535SBarry Smith 
8624b9ad928SBarry Smith /*@
86397177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
8644b9ad928SBarry Smith 
8654b9ad928SBarry Smith    Not Collective
8664b9ad928SBarry Smith 
8674b9ad928SBarry Smith    Input Parameter:
8684b9ad928SBarry Smith .  pc - the preconditioner context
8694b9ad928SBarry Smith 
8704b9ad928SBarry Smith    Output parameter:
8714b9ad928SBarry Smith .  levels - the number of levels
8724b9ad928SBarry Smith 
8734b9ad928SBarry Smith    Level: advanced
8744b9ad928SBarry Smith 
8754b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
8764b9ad928SBarry Smith 
87797177400SBarry Smith .seealso: PCMGSetLevels()
8784b9ad928SBarry Smith @*/
8797087cfbeSBarry Smith PetscErrorCode  PCMGGetLevels(PC pc,PetscInt *levels)
8804b9ad928SBarry Smith {
881f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
8824b9ad928SBarry Smith 
8834b9ad928SBarry Smith   PetscFunctionBegin;
8840700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
8854482741eSBarry Smith   PetscValidIntPointer(levels,2);
886f3fbd535SBarry Smith   *levels = mg->nlevels;
8874b9ad928SBarry Smith   PetscFunctionReturn(0);
8884b9ad928SBarry Smith }
8894b9ad928SBarry Smith 
8904b9ad928SBarry Smith /*@
89197177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
8924b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
8934b9ad928SBarry Smith 
894ad4df100SBarry Smith    Logically Collective on PC
8954b9ad928SBarry Smith 
8964b9ad928SBarry Smith    Input Parameters:
8974b9ad928SBarry Smith +  pc - the preconditioner context
8989dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
8999dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
9004b9ad928SBarry Smith 
9014b9ad928SBarry Smith    Options Database Key:
9024b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
9034b9ad928SBarry Smith    additive, full, kaskade
9044b9ad928SBarry Smith 
9054b9ad928SBarry Smith    Level: advanced
9064b9ad928SBarry Smith 
9074b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
9084b9ad928SBarry Smith 
90997177400SBarry Smith .seealso: PCMGSetLevels()
9104b9ad928SBarry Smith @*/
9117087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
9124b9ad928SBarry Smith {
913f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9144b9ad928SBarry Smith 
9154b9ad928SBarry Smith   PetscFunctionBegin;
9160700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
917c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
91831567311SBarry Smith   mg->am = form;
9199dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
920c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
921c60c7ad4SBarry Smith   PetscFunctionReturn(0);
922c60c7ad4SBarry Smith }
923c60c7ad4SBarry Smith 
924c60c7ad4SBarry Smith /*@
925c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
926c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
927c60c7ad4SBarry Smith 
928c60c7ad4SBarry Smith    Logically Collective on PC
929c60c7ad4SBarry Smith 
930c60c7ad4SBarry Smith    Input Parameter:
931c60c7ad4SBarry Smith .  pc - the preconditioner context
932c60c7ad4SBarry Smith 
933c60c7ad4SBarry Smith    Output Parameter:
934c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
935c60c7ad4SBarry Smith 
936c60c7ad4SBarry Smith 
937c60c7ad4SBarry Smith    Level: advanced
938c60c7ad4SBarry Smith 
939c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
940c60c7ad4SBarry Smith 
941c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
942c60c7ad4SBarry Smith @*/
943c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
944c60c7ad4SBarry Smith {
945c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
946c60c7ad4SBarry Smith 
947c60c7ad4SBarry Smith   PetscFunctionBegin;
948c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
949c60c7ad4SBarry Smith   *type = mg->am;
9504b9ad928SBarry Smith   PetscFunctionReturn(0);
9514b9ad928SBarry Smith }
9524b9ad928SBarry Smith 
9534b9ad928SBarry Smith /*@
9540d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
9554b9ad928SBarry Smith    complicated cycling.
9564b9ad928SBarry Smith 
957ad4df100SBarry Smith    Logically Collective on PC
9584b9ad928SBarry Smith 
9594b9ad928SBarry Smith    Input Parameters:
960c2be2410SBarry Smith +  pc - the multigrid context
961c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
9624b9ad928SBarry Smith 
9634b9ad928SBarry Smith    Options Database Key:
964c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
9654b9ad928SBarry Smith 
9664b9ad928SBarry Smith    Level: advanced
9674b9ad928SBarry Smith 
9684b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
9694b9ad928SBarry Smith 
9700d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
9714b9ad928SBarry Smith @*/
9727087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
9734b9ad928SBarry Smith {
974f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
975f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
97679416396SBarry Smith   PetscInt     i,levels;
9774b9ad928SBarry Smith 
9784b9ad928SBarry Smith   PetscFunctionBegin;
9790700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
98069fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
9811a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
982f3fbd535SBarry Smith   levels = mglevels[0]->levels;
9832fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
9844b9ad928SBarry Smith   PetscFunctionReturn(0);
9854b9ad928SBarry Smith }
9864b9ad928SBarry Smith 
9878cc2d5dfSBarry Smith /*@
9888cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
9898cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
9908cc2d5dfSBarry Smith 
991ad4df100SBarry Smith    Logically Collective on PC
9928cc2d5dfSBarry Smith 
9938cc2d5dfSBarry Smith    Input Parameters:
9948cc2d5dfSBarry Smith +  pc - the multigrid context
9958cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
9968cc2d5dfSBarry Smith 
9978cc2d5dfSBarry Smith    Options Database Key:
998e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
9998cc2d5dfSBarry Smith 
10008cc2d5dfSBarry Smith    Level: advanced
10018cc2d5dfSBarry Smith 
10028cc2d5dfSBarry Smith    Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
10038cc2d5dfSBarry Smith 
10048cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10058cc2d5dfSBarry Smith 
10068cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
10078cc2d5dfSBarry Smith @*/
10087087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
10098cc2d5dfSBarry Smith {
1010f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
10118cc2d5dfSBarry Smith 
10128cc2d5dfSBarry Smith   PetscFunctionBegin;
10130700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1014c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
10152a21e185SBarry Smith   mg->cyclesperpcapply = n;
10168cc2d5dfSBarry Smith   PetscFunctionReturn(0);
10178cc2d5dfSBarry Smith }
10188cc2d5dfSBarry Smith 
10192134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1020fb15c04eSBarry Smith {
1021fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1022fb15c04eSBarry Smith 
1023fb15c04eSBarry Smith   PetscFunctionBegin;
10242134b1e4SBarry Smith   mg->galerkin = use;
1025fb15c04eSBarry Smith   PetscFunctionReturn(0);
1026fb15c04eSBarry Smith }
1027fb15c04eSBarry Smith 
1028c2be2410SBarry Smith /*@
102997177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
103082b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1031c2be2410SBarry Smith 
1032ad4df100SBarry Smith    Logically Collective on PC
1033c2be2410SBarry Smith 
1034c2be2410SBarry Smith    Input Parameters:
1035c91913d3SJed Brown +  pc - the multigrid context
10362134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1037c2be2410SBarry Smith 
1038c2be2410SBarry Smith    Options Database Key:
10392134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1040c2be2410SBarry Smith 
1041c2be2410SBarry Smith    Level: intermediate
1042c2be2410SBarry Smith 
10431c1aac46SBarry Smith    Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
10441c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
10451c1aac46SBarry Smith 
1046c2be2410SBarry Smith .keywords: MG, set, Galerkin
1047c2be2410SBarry Smith 
10482134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
10493fc8bf9cSBarry Smith 
1050c2be2410SBarry Smith @*/
10512134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1052c2be2410SBarry Smith {
1053fb15c04eSBarry Smith   PetscErrorCode ierr;
1054c2be2410SBarry Smith 
1055c2be2410SBarry Smith   PetscFunctionBegin;
10560700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10572134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1058c2be2410SBarry Smith   PetscFunctionReturn(0);
1059c2be2410SBarry Smith }
1060c2be2410SBarry Smith 
10613fc8bf9cSBarry Smith /*@
10623fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
106382b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
10643fc8bf9cSBarry Smith 
10653fc8bf9cSBarry Smith    Not Collective
10663fc8bf9cSBarry Smith 
10673fc8bf9cSBarry Smith    Input Parameter:
10683fc8bf9cSBarry Smith .  pc - the multigrid context
10693fc8bf9cSBarry Smith 
10703fc8bf9cSBarry Smith    Output Parameter:
10712134b1e4SBarry 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
10723fc8bf9cSBarry Smith 
10733fc8bf9cSBarry Smith    Level: intermediate
10743fc8bf9cSBarry Smith 
10753fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
10763fc8bf9cSBarry Smith 
10772134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
10783fc8bf9cSBarry Smith 
10793fc8bf9cSBarry Smith @*/
10802134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
10813fc8bf9cSBarry Smith {
1082f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
10833fc8bf9cSBarry Smith 
10843fc8bf9cSBarry Smith   PetscFunctionBegin;
10850700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
10862134b1e4SBarry Smith   *galerkin = mg->galerkin;
10873fc8bf9cSBarry Smith   PetscFunctionReturn(0);
10883fc8bf9cSBarry Smith }
10893fc8bf9cSBarry Smith 
10904b9ad928SBarry Smith /*@
109106792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1092*710315b6SLawrence Mitchell    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1093*710315b6SLawrence Mitchell    pre- and post-smoothing steps.
109406792cafSBarry Smith 
109506792cafSBarry Smith    Logically Collective on PC
109606792cafSBarry Smith 
109706792cafSBarry Smith    Input Parameters:
109806792cafSBarry Smith +  mg - the multigrid context
109906792cafSBarry Smith -  n - the number of smoothing steps
110006792cafSBarry Smith 
110106792cafSBarry Smith    Options Database Key:
110206792cafSBarry Smith +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
110306792cafSBarry Smith 
110406792cafSBarry Smith    Level: advanced
110506792cafSBarry Smith 
110606792cafSBarry Smith    Notes: this does not set a value on the coarsest grid, since we assume that
110706792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
110806792cafSBarry Smith 
110906792cafSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
111006792cafSBarry Smith 
1111*710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp()
111206792cafSBarry Smith @*/
111306792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
111406792cafSBarry Smith {
111506792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
111606792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
111706792cafSBarry Smith   PetscErrorCode ierr;
111806792cafSBarry Smith   PetscInt       i,levels;
111906792cafSBarry Smith 
112006792cafSBarry Smith   PetscFunctionBegin;
112106792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
112206792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
11231a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
112406792cafSBarry Smith   levels = mglevels[0]->levels;
112506792cafSBarry Smith 
112606792cafSBarry Smith   for (i=1; i<levels; i++) {
112706792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
112806792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
112906792cafSBarry Smith     mg->default_smoothu = n;
113006792cafSBarry Smith     mg->default_smoothd = n;
113106792cafSBarry Smith   }
113206792cafSBarry Smith   PetscFunctionReturn(0);
113306792cafSBarry Smith }
113406792cafSBarry Smith 
1135f442ab6aSBarry Smith /*@
1136f442ab6aSBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels
1137*710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1138f442ab6aSBarry Smith 
1139f442ab6aSBarry Smith    Logically Collective on PC
1140f442ab6aSBarry Smith 
1141f442ab6aSBarry Smith    Input Parameters:
1142f442ab6aSBarry Smith .  pc - the preconditioner context
1143f442ab6aSBarry Smith 
1144f442ab6aSBarry Smith    Options Database Key:
1145f442ab6aSBarry Smith .  -pc_mg_distinct_smoothup
1146f442ab6aSBarry Smith 
1147f442ab6aSBarry Smith    Level: advanced
1148f442ab6aSBarry Smith 
1149f442ab6aSBarry Smith    Notes: this does not set a value on the coarsest grid, since we assume that
1150f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1151f442ab6aSBarry Smith 
1152f442ab6aSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1153f442ab6aSBarry Smith 
1154*710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth()
1155f442ab6aSBarry Smith @*/
1156f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1157f442ab6aSBarry Smith {
1158f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1159f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1160f442ab6aSBarry Smith   PetscErrorCode ierr;
1161f442ab6aSBarry Smith   PetscInt       i,levels;
1162f442ab6aSBarry Smith   KSP            subksp;
1163f442ab6aSBarry Smith 
1164f442ab6aSBarry Smith   PetscFunctionBegin;
1165f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1166f442ab6aSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1167f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1168f442ab6aSBarry Smith 
1169f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1170*710315b6SLawrence Mitchell     const char *prefix = NULL;
1171f442ab6aSBarry Smith     /* make sure smoother up and down are different */
1172f442ab6aSBarry Smith     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1173*710315b6SLawrence Mitchell     ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr);
1174*710315b6SLawrence Mitchell     ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr);
1175f442ab6aSBarry Smith     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1176f442ab6aSBarry Smith   }
1177f442ab6aSBarry Smith   PetscFunctionReturn(0);
1178f442ab6aSBarry Smith }
1179f442ab6aSBarry Smith 
11804b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
11814b9ad928SBarry Smith 
11823b09bd56SBarry Smith /*MC
1183ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
11843b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
11853b09bd56SBarry Smith 
11863b09bd56SBarry Smith    Options Database Keys:
11873b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
11884afba7b0SPatrick Sanan .  -pc_mg_cycle_type <v,w> -
11898c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
11903b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1191*710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
11922134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
11938cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
11948cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1195e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
11968cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
11978cf6037eSBarry Smith                         to the binary output file called binaryoutput
11983b09bd56SBarry Smith 
1199e941fc33SBarry 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
12003b09bd56SBarry Smith 
12018cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
12028cf6037eSBarry Smith 
120323067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
120423067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
120523067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
120623067569SBarry Smith        residual is computed at the end of each cycle.
120723067569SBarry Smith 
12083b09bd56SBarry Smith    Level: intermediate
12093b09bd56SBarry Smith 
12108f87f92bSBarry Smith    Concepts: multigrid/multilevel
12113b09bd56SBarry Smith 
12128cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1213*710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1214*710315b6SLawrence Mitchell            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
121597177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12160d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12173b09bd56SBarry Smith M*/
12183b09bd56SBarry Smith 
12198cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
12204b9ad928SBarry Smith {
1221f3fbd535SBarry Smith   PC_MG          *mg;
1222f3fbd535SBarry Smith   PetscErrorCode ierr;
1223f3fbd535SBarry Smith 
12244b9ad928SBarry Smith   PetscFunctionBegin;
1225b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1226f3fbd535SBarry Smith   pc->data     = (void*)mg;
1227f3fbd535SBarry Smith   mg->nlevels  = -1;
122810eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
12292134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1230f3fbd535SBarry Smith 
123137a44384SMark Adams   pc->useAmat = PETSC_TRUE;
123237a44384SMark Adams 
12334b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
12344b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1235a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
12364b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
12374b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
12384b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1239fb15c04eSBarry Smith 
1240fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
12414b9ad928SBarry Smith   PetscFunctionReturn(0);
12424b9ad928SBarry Smith }
1243