xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision 37b5128c94a41d620f47f34aadb272941d4da11e)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith     Defines the multigrid preconditioner interface.
44b9ad928SBarry Smith */
5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>                    /*I "petscksp.h" I*/
61e25c274SJed Brown #include <petscdm.h>
7391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*);
84b9ad928SBarry Smith 
931567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
104b9ad928SBarry Smith {
1131567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1231567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
136849ba73SBarry Smith   PetscErrorCode ierr;
1431567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
15a0808db4SHong Zhang   PC             subpc;
16a0808db4SHong Zhang   PCFailedReason pcreason;
174b9ad928SBarry Smith 
184b9ad928SBarry Smith   PetscFunctionBegin;
1963e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2031567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
21a0808db4SHong Zhang   ierr = KSPGetPC(mglevels->smoothd,&subpc);CHKERRQ(ierr);
22a0808db4SHong Zhang   ierr = PCGetSetUpFailedReason(subpc,&pcreason);CHKERRQ(ierr);
23a0808db4SHong Zhang   if (pcreason) {
24a0808db4SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
25a0808db4SHong Zhang   }
2663e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2731567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2863e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2931567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
3063e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
314b9ad928SBarry Smith 
324b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
3331567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
344b9ad928SBarry Smith       PetscReal rnorm;
3531567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
364b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3770441072SBarry Smith         if (rnorm < mg->abstol) {
384d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
3957622a8eSBarry 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);
404b9ad928SBarry Smith         } else {
414d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
4257622a8eSBarry 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);
434b9ad928SBarry Smith         }
444b9ad928SBarry Smith         PetscFunctionReturn(0);
454b9ad928SBarry Smith       }
464b9ad928SBarry Smith     }
474b9ad928SBarry Smith 
4831567311SBarry Smith     mgc = *(mglevelsin - 1);
4963e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5031567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
5163e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
52efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
534b9ad928SBarry Smith     while (cycles--) {
5431567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
554b9ad928SBarry Smith     }
5663e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5731567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5863e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5963e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
6031567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
6163e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
624b9ad928SBarry Smith   }
634b9ad928SBarry Smith   PetscFunctionReturn(0);
644b9ad928SBarry Smith }
654b9ad928SBarry Smith 
66ace3abfcSBarry 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)
674b9ad928SBarry Smith {
68f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
69f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
70dfbe8321SBarry Smith   PetscErrorCode ierr;
71391689abSStefano Zampini   PC             tpc;
72391689abSStefano Zampini   PetscBool      changeu,changed;
73f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
744b9ad928SBarry Smith 
754b9ad928SBarry Smith   PetscFunctionBegin;
76a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
77a762d673SBarry Smith   for (i=0; i<levels; i++) {
78a762d673SBarry Smith     if (!mglevels[i]->A) {
79a762d673SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
80a762d673SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
81a762d673SBarry Smith     }
82a762d673SBarry Smith   }
83391689abSStefano Zampini 
84391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
85391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
86391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
87391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
88391689abSStefano Zampini   if (!changed && !changeu) {
89391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
90f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
91391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
92391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
93391689abSStefano Zampini       Vec *vec;
94391689abSStefano Zampini 
95391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
96391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
977ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
98391689abSStefano Zampini     }
99391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
100391689abSStefano Zampini   }
101f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
1024b9ad928SBarry Smith 
10331567311SBarry Smith   mg->rtol   = rtol;
10431567311SBarry Smith   mg->abstol = abstol;
10531567311SBarry Smith   mg->dtol   = dtol;
1064b9ad928SBarry Smith   if (rtol) {
1074b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
1084b9ad928SBarry Smith     PetscReal rnorm;
1097319c654SBarry Smith     if (zeroguess) {
1107319c654SBarry Smith       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
1117319c654SBarry Smith     } else {
112f3fbd535SBarry Smith       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
1134b9ad928SBarry Smith       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
1147319c654SBarry Smith     }
11531567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
1162fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
1172fa5cd67SKarl Rupp   else mg->ttol = 0.0;
1184b9ad928SBarry Smith 
1194d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
120336babb1SJed Brown      stop prematurely due to small residual */
1214d0a8057SBarry Smith   for (i=1; i<levels; i++) {
122f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
123f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
12423067569SBarry Smith       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
12523067569SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
126f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1274b9ad928SBarry Smith     }
1284d0a8057SBarry Smith   }
1294d0a8057SBarry Smith 
1304d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1314d0a8057SBarry Smith   for (i=0; i<its; i++) {
132f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1334d0a8057SBarry Smith     if (*reason) break;
1344d0a8057SBarry Smith   }
1354d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1364d0a8057SBarry Smith   *outits = i;
137391689abSStefano Zampini   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
1384b9ad928SBarry Smith   PetscFunctionReturn(0);
1394b9ad928SBarry Smith }
1404b9ad928SBarry Smith 
1413aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1423aeaf226SBarry Smith {
1433aeaf226SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1443aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1453aeaf226SBarry Smith   PetscErrorCode ierr;
1463aeaf226SBarry Smith   PetscInt       i,n;
1473aeaf226SBarry Smith 
1483aeaf226SBarry Smith   PetscFunctionBegin;
1493aeaf226SBarry Smith   if (mglevels) {
1503aeaf226SBarry Smith     n = mglevels[0]->levels;
1513aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1523aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1533aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1543aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
1553aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1563aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
1570de8493bSLawrence Mitchell       ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr);
15873250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1593aeaf226SBarry Smith     }
160391689abSStefano Zampini     /* this is not null only if the smoother on the finest level
161391689abSStefano Zampini        changes the rhs during PreSolve */
162391689abSStefano Zampini     ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr);
1633aeaf226SBarry Smith 
1643aeaf226SBarry Smith     for (i=0; i<n; i++) {
1653aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1663aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1673aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1683aeaf226SBarry Smith       }
1693aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1703aeaf226SBarry Smith     }
1713aeaf226SBarry Smith   }
1723aeaf226SBarry Smith   PetscFunctionReturn(0);
1733aeaf226SBarry Smith }
1743aeaf226SBarry Smith 
17597d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
1764b9ad928SBarry Smith {
177dfbe8321SBarry Smith   PetscErrorCode ierr;
178f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
179ce94432eSBarry Smith   MPI_Comm       comm;
1803aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
18110eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
18210167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
183f3fbd535SBarry Smith   PetscInt       i;
184f3fbd535SBarry Smith   PetscMPIInt    size;
185f3fbd535SBarry Smith   const char     *prefix;
186f3fbd535SBarry Smith   PC             ipc;
1873aeaf226SBarry Smith   PetscInt       n;
1884b9ad928SBarry Smith 
1894b9ad928SBarry Smith   PetscFunctionBegin;
1900700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
191548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
192548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1931a97d7b6SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1943aeaf226SBarry Smith   if (mglevels) {
19510eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
1963aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
1973aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
1983aeaf226SBarry Smith     n    = mglevels[0]->levels;
1993aeaf226SBarry Smith     for (i=0; i<n; i++) {
2003aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2013aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
2023aeaf226SBarry Smith       }
2033aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
2043aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
2053aeaf226SBarry Smith     }
2063aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
2073aeaf226SBarry Smith   }
208f3fbd535SBarry Smith 
209f3fbd535SBarry Smith   mg->nlevels = levels;
210f3fbd535SBarry Smith 
211785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
2123bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
213f3fbd535SBarry Smith 
214f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
215f3fbd535SBarry Smith 
216a9db87a2SMatthew G Knepley   mg->stageApply = 0;
217f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
218b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
2192fa5cd67SKarl Rupp 
220f3fbd535SBarry Smith     mglevels[i]->level               = i;
221f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
22210eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
223336babb1SJed Brown     mg->default_smoothu              = 2;
224336babb1SJed Brown     mg->default_smoothd              = 2;
22563e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
22663e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
22763e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
22863e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
229f3fbd535SBarry Smith 
230f3fbd535SBarry Smith     if (comms) comm = comms[i];
231f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
232422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
2330ee9a668SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
2340ee9a668SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
2350ee9a668SBarry Smith     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
2360ee9a668SBarry Smith     if (i || levels == 1) {
2370ee9a668SBarry Smith       char tprefix[128];
2380ee9a668SBarry Smith 
239336babb1SJed Brown       ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2400059c7bdSJed Brown       ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
241336babb1SJed Brown       ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
242336babb1SJed Brown       ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
243336babb1SJed Brown       ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
2440ee9a668SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
245f3fbd535SBarry Smith 
2460ee9a668SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
2470ee9a668SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
2480ee9a668SBarry Smith     } else {
249f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
250f3fbd535SBarry Smith 
2519dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
252f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
253f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
254f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
255f3fbd535SBarry Smith       if (size > 1) {
256f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
257f3fbd535SBarry Smith       } else {
258f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
259f3fbd535SBarry Smith       }
260753b7fb9SBarry Smith       ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
261f3fbd535SBarry Smith     }
2623bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
2632fa5cd67SKarl Rupp 
264f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
26531567311SBarry Smith     mg->rtol             = 0.0;
26631567311SBarry Smith     mg->abstol           = 0.0;
26731567311SBarry Smith     mg->dtol             = 0.0;
26831567311SBarry Smith     mg->ttol             = 0.0;
26931567311SBarry Smith     mg->cyclesperpcapply = 1;
270f3fbd535SBarry Smith   }
271f3fbd535SBarry Smith   mg->levels = mglevels;
27210eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2734b9ad928SBarry Smith   PetscFunctionReturn(0);
2744b9ad928SBarry Smith }
2754b9ad928SBarry Smith 
27697d33e41SMatthew G. Knepley /*@C
27797d33e41SMatthew G. Knepley    PCMGSetLevels - Sets the number of levels to use with MG.
27897d33e41SMatthew G. Knepley    Must be called before any other MG routine.
27997d33e41SMatthew G. Knepley 
28097d33e41SMatthew G. Knepley    Logically Collective on PC
28197d33e41SMatthew G. Knepley 
28297d33e41SMatthew G. Knepley    Input Parameters:
28397d33e41SMatthew G. Knepley +  pc - the preconditioner context
28497d33e41SMatthew G. Knepley .  levels - the number of levels
28597d33e41SMatthew G. Knepley -  comms - optional communicators for each level; this is to allow solving the coarser problems
28697d33e41SMatthew G. Knepley            on smaller sets of processors.
28797d33e41SMatthew G. Knepley 
28897d33e41SMatthew G. Knepley    Level: intermediate
28997d33e41SMatthew G. Knepley 
29097d33e41SMatthew G. Knepley    Notes:
29197d33e41SMatthew G. Knepley      If the number of levels is one then the multigrid uses the -mg_levels prefix
29297d33e41SMatthew G. Knepley   for setting the level options rather than the -mg_coarse prefix.
29397d33e41SMatthew G. Knepley 
29497d33e41SMatthew G. Knepley .keywords: MG, set, levels, multigrid
29597d33e41SMatthew G. Knepley 
29697d33e41SMatthew G. Knepley .seealso: PCMGSetType(), PCMGGetLevels()
29797d33e41SMatthew G. Knepley @*/
29897d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
29997d33e41SMatthew G. Knepley {
30097d33e41SMatthew G. Knepley   PetscErrorCode ierr;
30197d33e41SMatthew G. Knepley 
30297d33e41SMatthew G. Knepley   PetscFunctionBegin;
30397d33e41SMatthew G. Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
30497d33e41SMatthew G. Knepley   if (comms) PetscValidPointer(comms,3);
305*37b5128cSMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr);
30697d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
30797d33e41SMatthew G. Knepley }
30897d33e41SMatthew G. Knepley 
309c07bf074SBarry Smith 
310c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
311c07bf074SBarry Smith {
312c07bf074SBarry Smith   PetscErrorCode ierr;
313a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
314a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
315a06653b4SBarry Smith   PetscInt       i,n;
316c07bf074SBarry Smith 
317c07bf074SBarry Smith   PetscFunctionBegin;
318a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
319a06653b4SBarry Smith   if (mglevels) {
320a06653b4SBarry Smith     n = mglevels[0]->levels;
321a06653b4SBarry Smith     for (i=0; i<n; i++) {
322a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
3236bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
324a06653b4SBarry Smith       }
3256bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
326a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
327a06653b4SBarry Smith     }
328a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
329a06653b4SBarry Smith   }
330c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
331f3fbd535SBarry Smith   PetscFunctionReturn(0);
332f3fbd535SBarry Smith }
333f3fbd535SBarry Smith 
334f3fbd535SBarry Smith 
335f3fbd535SBarry Smith 
33609573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
33709573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
33809573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
339f3fbd535SBarry Smith 
340f3fbd535SBarry Smith /*
341f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
342f3fbd535SBarry Smith              or full cycle of multigrid.
343f3fbd535SBarry Smith 
344f3fbd535SBarry Smith   Note:
345f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
346f3fbd535SBarry Smith */
347f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
348f3fbd535SBarry Smith {
349f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
350f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
351f3fbd535SBarry Smith   PetscErrorCode ierr;
352391689abSStefano Zampini   PC             tpc;
353f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
354391689abSStefano Zampini   PetscBool      changeu,changed;
355f3fbd535SBarry Smith 
356f3fbd535SBarry Smith   PetscFunctionBegin;
357a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
358e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
359298cc208SBarry Smith   for (i=0; i<levels; i++) {
360e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
36123ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
362298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
363e1d8e5deSBarry Smith     }
364e1d8e5deSBarry Smith   }
365e1d8e5deSBarry Smith 
366391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
367391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
368391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
369391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
370391689abSStefano Zampini   if (!changeu && !changed) {
371391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
372f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
373391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
374391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
375391689abSStefano Zampini       Vec *vec;
376391689abSStefano Zampini 
377391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
378391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
3797ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
380391689abSStefano Zampini     }
381391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
382391689abSStefano Zampini   }
383f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
384391689abSStefano Zampini 
38531567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
386f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
38731567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
3880298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
389f3fbd535SBarry Smith     }
3902fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
39131567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
3922fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
39331567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
3942fa5cd67SKarl Rupp   } else {
395f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
396f3fbd535SBarry Smith   }
397a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
398391689abSStefano Zampini   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
399f3fbd535SBarry Smith   PetscFunctionReturn(0);
400f3fbd535SBarry Smith }
401f3fbd535SBarry Smith 
402f3fbd535SBarry Smith 
4034416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
404f3fbd535SBarry Smith {
405f3fbd535SBarry Smith   PetscErrorCode   ierr;
406710315b6SLawrence Mitchell   PetscInt         levels,cycles;
4072134b1e4SBarry Smith   PetscBool        flg;
408f3fbd535SBarry Smith   PC_MG            *mg = (PC_MG*)pc->data;
4093d3eaba7SBarry Smith   PC_MG_Levels     **mglevels;
410f3fbd535SBarry Smith   PCMGType         mgtype;
411f3fbd535SBarry Smith   PCMGCycleType    mgctype;
4122134b1e4SBarry Smith   PCMGGalerkinType gtype;
413f3fbd535SBarry Smith 
414f3fbd535SBarry Smith   PetscFunctionBegin;
4151d743356SStefano Zampini   levels = PetscMax(mg->nlevels,1);
416e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
417f3fbd535SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
4181a97d7b6SStefano Zampini   if (!flg && !mg->levels && pc->dm) {
419eb3f98d2SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
420eb3f98d2SBarry Smith     levels++;
4213aeaf226SBarry Smith     mg->usedmfornumberoflevels = PETSC_TRUE;
422eb3f98d2SBarry Smith   }
4230298fd71SBarry Smith   ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
4243aeaf226SBarry Smith   mglevels = mg->levels;
4253aeaf226SBarry Smith 
426f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
427f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
428f3fbd535SBarry Smith   if (flg) {
429f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
4302fa5cd67SKarl Rupp   }
4312134b1e4SBarry Smith   gtype = mg->galerkin;
4322134b1e4SBarry Smith   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
4332134b1e4SBarry Smith   if (flg) {
4342134b1e4SBarry Smith     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
435f3fbd535SBarry Smith   }
436f56b1265SBarry Smith   flg = PETSC_FALSE;
437f442ab6aSBarry Smith   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
438f442ab6aSBarry Smith   if (flg) {
439f442ab6aSBarry Smith     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
440f442ab6aSBarry Smith   }
44131567311SBarry Smith   mgtype = mg->am;
442f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
443f3fbd535SBarry Smith   if (flg) {
444f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
445f3fbd535SBarry Smith   }
44631567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
4475363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
448f3fbd535SBarry Smith     if (flg) {
449f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
450f3fbd535SBarry Smith     }
451f3fbd535SBarry Smith   }
452f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4530298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
454f3fbd535SBarry Smith   if (flg) {
455f3fbd535SBarry Smith     PetscInt i;
456f3fbd535SBarry Smith     char     eventname[128];
4571a97d7b6SStefano Zampini 
458f3fbd535SBarry Smith     levels = mglevels[0]->levels;
459f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
460f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
46163e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
462f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
46363e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
464f3fbd535SBarry Smith       if (i) {
465f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
46663e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
467f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
46863e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
469f3fbd535SBarry Smith       }
470f3fbd535SBarry Smith     }
471eec35531SMatthew G Knepley 
472ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
473eec35531SMatthew G Knepley     {
474eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
475eec35531SMatthew G Knepley       PetscStageLog stageLog;
476eec35531SMatthew G Knepley       PetscInt      st;
477eec35531SMatthew G Knepley 
478eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
479eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
480eec35531SMatthew G Knepley         PetscBool same;
481eec35531SMatthew G Knepley 
482eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4832fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
484eec35531SMatthew G Knepley       }
485eec35531SMatthew G Knepley       if (!mg->stageApply) {
486eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
487eec35531SMatthew G Knepley       }
488eec35531SMatthew G Knepley     }
489ec60ca73SSatish Balay #endif
490f3fbd535SBarry Smith   }
491f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
492f3fbd535SBarry Smith   PetscFunctionReturn(0);
493f3fbd535SBarry Smith }
494f3fbd535SBarry Smith 
4956a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4966a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
4972134b1e4SBarry Smith const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0};
498f3fbd535SBarry Smith 
4999804daf3SBarry Smith #include <petscdraw.h>
500f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
501f3fbd535SBarry Smith {
502f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
503f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
504f3fbd535SBarry Smith   PetscErrorCode ierr;
505e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
506219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
507f3fbd535SBarry Smith 
508f3fbd535SBarry Smith   PetscFunctionBegin;
509251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5105b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
511219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
512f3fbd535SBarry Smith   if (iascii) {
513e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
514efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
51531567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
51631567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
517f3fbd535SBarry Smith     }
5182134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
519f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
5202134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
5212134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
5222134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
5232134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
5242134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
5252134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
5264f66f45eSBarry Smith     } else {
5274f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
528f3fbd535SBarry Smith     }
5295adeb434SBarry Smith     if (mg->view){
5305adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
5315adeb434SBarry Smith     }
532f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
533f3fbd535SBarry Smith       if (!i) {
53489cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
535f3fbd535SBarry Smith       } else {
53689cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
537f3fbd535SBarry Smith       }
538f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
539f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
540f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
541f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
542f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
543f3fbd535SBarry Smith       } else if (i) {
54489cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
545f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
546f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
547f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
548f3fbd535SBarry Smith       }
549f3fbd535SBarry Smith     }
5505b0b0462SBarry Smith   } else if (isbinary) {
5515b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5525b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5535b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5545b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5555b0b0462SBarry Smith       }
5565b0b0462SBarry Smith     }
557219b1045SBarry Smith   } else if (isdraw) {
558219b1045SBarry Smith     PetscDraw draw;
559b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
560219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
561219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5620298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
563219b1045SBarry Smith     bottom = y - th;
564219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
565b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
566219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
567219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
568219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
569b4375e8dSPeter Brune       } else {
570b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
571b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
572b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
573b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
574b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
575b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
576b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
577b4375e8dSPeter Brune       }
5780298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5791cd381d6SBarry Smith       bottom -= th;
580219b1045SBarry Smith     }
5815b0b0462SBarry Smith   }
582f3fbd535SBarry Smith   PetscFunctionReturn(0);
583f3fbd535SBarry Smith }
584f3fbd535SBarry Smith 
585af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
586af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
587cab2e9ccSBarry Smith 
588f3fbd535SBarry Smith /*
589f3fbd535SBarry Smith     Calls setup for the KSP on each level
590f3fbd535SBarry Smith */
591f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
592f3fbd535SBarry Smith {
593f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
594f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
595f3fbd535SBarry Smith   PetscErrorCode ierr;
5961a97d7b6SStefano Zampini   PetscInt       i,n;
59798e04984SBarry Smith   PC             cpc;
5983db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
599f3fbd535SBarry Smith   Mat            dA,dB;
600f3fbd535SBarry Smith   Vec            tvec;
601218a07d4SBarry Smith   DM             *dms;
602649052a6SBarry Smith   PetscViewer    viewer = 0;
6032134b1e4SBarry Smith   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
604f3fbd535SBarry Smith 
605f3fbd535SBarry Smith   PetscFunctionBegin;
6061a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
6071a97d7b6SStefano Zampini   n = mglevels[0]->levels;
60801bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
6093aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
6103aeaf226SBarry Smith     PetscInt levels;
6113aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
6123aeaf226SBarry Smith     levels++;
6133aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
6140298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
6153aeaf226SBarry Smith       n        = levels;
6163aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
6173aeaf226SBarry Smith       mglevels =  mg->levels;
6183aeaf226SBarry Smith     }
6193aeaf226SBarry Smith   }
62098e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
6213aeaf226SBarry Smith 
622f3fbd535SBarry Smith 
623f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
624f3fbd535SBarry Smith   /* so use those from global PC */
625f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
6260298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
62798e04984SBarry Smith   if (opsset) {
62898e04984SBarry Smith     Mat mmat;
62923ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
63098e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
63198e04984SBarry Smith   }
6324a5f07a5SMark F. Adams 
6334a5f07a5SMark F. Adams   if (!opsset) {
63471b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
6354a5f07a5SMark F. Adams     if(use_amat){
636f3fbd535SBarry 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);
63723ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
638f3fbd535SBarry Smith     }
6394a5f07a5SMark F. Adams     else {
6404a5f07a5SMark 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);
64123ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
6424a5f07a5SMark F. Adams     }
6434a5f07a5SMark F. Adams   }
644f3fbd535SBarry Smith 
6459c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
6469c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
6479c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
6489c7ffb39SBarry Smith       continue;
6499c7ffb39SBarry Smith     }
6509c7ffb39SBarry Smith   }
6512134b1e4SBarry Smith 
6522134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
6532134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
6542134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
6552134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
6562134b1e4SBarry Smith   }
6572134b1e4SBarry Smith 
6582134b1e4SBarry Smith 
6599c7ffb39SBarry Smith   /*
6609c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6619c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6629c7ffb39SBarry Smith   */
6632134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
6642d2b81a6SBarry Smith 	/* construct the interpolation from the DMs */
6652e499ae9SBarry Smith     Mat p;
66673250ac0SBarry Smith     Vec rscale;
667785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
668218a07d4SBarry Smith     dms[n-1] = pc->dm;
669ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
670ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
67141fce8fdSHong Zhang 	/*
67241fce8fdSHong Zhang 	   Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
67341fce8fdSHong Zhang 	   Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
67441fce8fdSHong Zhang 	   But it is safe to use -dm_mat_type.
67541fce8fdSHong Zhang 
67641fce8fdSHong Zhang 	   The mat type should not be hardcoded like this, we need to find a better way.
67741fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
67841fce8fdSHong Zhang     */
679218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
680942e3340SBarry Smith       DMKSP     kdm;
681eab52d2bSLawrence Mitchell       PetscBool dmhasrestrict, dmhasinject;
6823c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6832134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
684942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
685d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
686d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6870298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6880298fd71SBarry Smith       kdm->rhsctx          = NULL;
68924c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
690e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6912d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
69200ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
69373250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6946bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
695218a07d4SBarry Smith       }
6963ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
6973ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
6983ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
6993ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
7003ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
7013ad4599aSBarry Smith       }
702eab52d2bSLawrence Mitchell       ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr);
703eab52d2bSLawrence Mitchell       if (dmhasinject && !mglevels[i+1]->inject){
704eab52d2bSLawrence Mitchell         ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr);
705eab52d2bSLawrence Mitchell         ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr);
706eab52d2bSLawrence Mitchell         ierr = MatDestroy(&p);CHKERRQ(ierr);
707eab52d2bSLawrence Mitchell       }
70824c3aa18SBarry Smith     }
7092d2b81a6SBarry Smith 
710ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
7112d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
712ce4cda84SJed Brown   }
713cab2e9ccSBarry Smith 
714ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
7152134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
716cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
717cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
718218a07d4SBarry Smith   }
719218a07d4SBarry Smith 
7202134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
7212134b1e4SBarry Smith     Mat       A,B;
7222134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
7232134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
7242134b1e4SBarry Smith 
7252134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
7262134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
7272134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
728f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
729b9367914SBarry 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");
730b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
731b9367914SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
732b9367914SBarry Smith       }
733b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
734b9367914SBarry Smith         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
735b9367914SBarry Smith       }
7362134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
7372134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
7382134b1e4SBarry Smith       }
7392134b1e4SBarry Smith       if (doA) {
7402df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
7412134b1e4SBarry Smith       }
7422134b1e4SBarry Smith       if (doB) {
7432df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
744b9367914SBarry Smith       }
7452134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
7462134b1e4SBarry Smith       if (!doA && dAeqdB) {
7472134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
7482134b1e4SBarry Smith         A = B;
7492134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX ) {
7502134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
7512134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
752b9367914SBarry Smith       }
7532134b1e4SBarry Smith       if (!doB && dAeqdB) {
7542134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
7552134b1e4SBarry Smith         B = A;
7562134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
75723ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
7582134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
7597d537d92SJed Brown       }
7602134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
7612134b1e4SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
7622134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
7632134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
7642134b1e4SBarry Smith       }
7652134b1e4SBarry Smith       dA = A;
766f3fbd535SBarry Smith       dB = B;
767f3fbd535SBarry Smith     }
768f3fbd535SBarry Smith   }
7692134b1e4SBarry Smith   if (needRestricts && pc->dm && pc->dm->x) {
770cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
771cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
772c88c5224SJed Brown       Mat R;
773c88c5224SJed Brown       Vec rscale;
774cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
775cab2e9ccSBarry Smith         Vec *vecs;
7762a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
777cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
778cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
779cab2e9ccSBarry Smith       }
780c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
781c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
782c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
783c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
784cab2e9ccSBarry Smith     }
785f3fbd535SBarry Smith   }
7862134b1e4SBarry Smith   if (needRestricts && pc->dm) {
787caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
788ccceb50cSJed Brown       DM  dmfine,dmcoarse;
789caa4e7f2SJed Brown       Mat Restrict,Inject;
790caa4e7f2SJed Brown       Vec rscale;
791ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
792ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
793caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
794caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
795eab52d2bSLawrence Mitchell       ierr   = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr);
796caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
797caa4e7f2SJed Brown     }
798caa4e7f2SJed Brown   }
799f3fbd535SBarry Smith 
800f3fbd535SBarry Smith   if (!pc->setupcalled) {
801f3fbd535SBarry Smith     for (i=0; i<n; i++) {
802f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
803f3fbd535SBarry Smith     }
804f3fbd535SBarry Smith     for (i=1; i<n; i++) {
805f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
806f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
807f3fbd535SBarry Smith       }
808f3fbd535SBarry Smith     }
8093ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
810f3fbd535SBarry Smith     for (i=1; i<n; i++) {
8113ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
8123ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
813f3fbd535SBarry Smith     }
814f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
815f3fbd535SBarry Smith       if (!mglevels[i]->b) {
816f3fbd535SBarry Smith         Vec *vec;
8172a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
818f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
8196bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
820f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
821f3fbd535SBarry Smith       }
822f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
823f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
824f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
8256bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
826f3fbd535SBarry Smith       }
827f3fbd535SBarry Smith       if (!mglevels[i]->x) {
828f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
829f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
8306bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
831f3fbd535SBarry Smith       }
832f3fbd535SBarry Smith     }
833f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
834f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
835f3fbd535SBarry Smith       Vec *vec;
8362a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
837f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
8386bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
839f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
840f3fbd535SBarry Smith     }
841f3fbd535SBarry Smith   }
842f3fbd535SBarry Smith 
84398e04984SBarry Smith   if (pc->dm) {
84498e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
84598e04984SBarry Smith     for (i=0; i<n-1; i++) {
84698e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
84798e04984SBarry Smith     }
84898e04984SBarry Smith   }
849f3fbd535SBarry Smith 
850f3fbd535SBarry Smith   for (i=1; i<n; i++) {
8512a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
852f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
853f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
854f3fbd535SBarry Smith     }
85563e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
856f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
857899639b0SHong Zhang     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
858899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
859899639b0SHong Zhang     }
86063e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
861d42688cbSBarry Smith     if (!mglevels[i]->residual) {
862d42688cbSBarry Smith       Mat mat;
86313842ffbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
86454b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
865d42688cbSBarry Smith     }
866f3fbd535SBarry Smith   }
867f3fbd535SBarry Smith   for (i=1; i<n; i++) {
868f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
869f3fbd535SBarry Smith       Mat downmat,downpmat;
870f3fbd535SBarry Smith 
871f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
8720298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
873f3fbd535SBarry Smith       if (!opsset) {
87423ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
87523ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
876f3fbd535SBarry Smith       }
877f3fbd535SBarry Smith 
878f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
87963e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
880f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
881899639b0SHong Zhang       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) {
882899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
883899639b0SHong Zhang       }
88463e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
885f3fbd535SBarry Smith     }
886f3fbd535SBarry Smith   }
887f3fbd535SBarry Smith 
88863e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
889f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
890899639b0SHong Zhang   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) {
891899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
892899639b0SHong Zhang   }
89363e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
894f3fbd535SBarry Smith 
895f3fbd535SBarry Smith   /*
896f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
897e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
898f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
899f3fbd535SBarry Smith 
900f3fbd535SBarry Smith    Only support one or the other at the same time.
901f3fbd535SBarry Smith   */
902f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
903c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
904ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
905f3fbd535SBarry Smith   dump = PETSC_FALSE;
906f3fbd535SBarry Smith #endif
907c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
908ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
909f3fbd535SBarry Smith 
910f3fbd535SBarry Smith   if (viewer) {
911f3fbd535SBarry Smith     for (i=1; i<n; i++) {
912f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
913f3fbd535SBarry Smith     }
914f3fbd535SBarry Smith     for (i=0; i<n; i++) {
915f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
916f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
917f3fbd535SBarry Smith     }
918f3fbd535SBarry Smith   }
919f3fbd535SBarry Smith   PetscFunctionReturn(0);
920f3fbd535SBarry Smith }
921f3fbd535SBarry Smith 
922f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
923f3fbd535SBarry Smith 
92497d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
92597d33e41SMatthew G. Knepley {
92697d33e41SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
92797d33e41SMatthew G. Knepley 
92897d33e41SMatthew G. Knepley   PetscFunctionBegin;
92997d33e41SMatthew G. Knepley   *levels = mg->nlevels;
93097d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
93197d33e41SMatthew G. Knepley }
93297d33e41SMatthew G. Knepley 
9334b9ad928SBarry Smith /*@
93497177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
9354b9ad928SBarry Smith 
9364b9ad928SBarry Smith    Not Collective
9374b9ad928SBarry Smith 
9384b9ad928SBarry Smith    Input Parameter:
9394b9ad928SBarry Smith .  pc - the preconditioner context
9404b9ad928SBarry Smith 
9414b9ad928SBarry Smith    Output parameter:
9424b9ad928SBarry Smith .  levels - the number of levels
9434b9ad928SBarry Smith 
9444b9ad928SBarry Smith    Level: advanced
9454b9ad928SBarry Smith 
9464b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
9474b9ad928SBarry Smith 
94897177400SBarry Smith .seealso: PCMGSetLevels()
9494b9ad928SBarry Smith @*/
9507087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
9514b9ad928SBarry Smith {
952e952c529SMatthew G. Knepley   PetscErrorCode ierr;
9534b9ad928SBarry Smith 
9544b9ad928SBarry Smith   PetscFunctionBegin;
9550700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9564482741eSBarry Smith   PetscValidIntPointer(levels,2);
95797d33e41SMatthew G. Knepley   *levels = 0;
958*37b5128cSMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr);
9594b9ad928SBarry Smith   PetscFunctionReturn(0);
9604b9ad928SBarry Smith }
9614b9ad928SBarry Smith 
9624b9ad928SBarry Smith /*@
96397177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
9644b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
9654b9ad928SBarry Smith 
966ad4df100SBarry Smith    Logically Collective on PC
9674b9ad928SBarry Smith 
9684b9ad928SBarry Smith    Input Parameters:
9694b9ad928SBarry Smith +  pc - the preconditioner context
9709dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
9719dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
9724b9ad928SBarry Smith 
9734b9ad928SBarry Smith    Options Database Key:
9744b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
9754b9ad928SBarry Smith    additive, full, kaskade
9764b9ad928SBarry Smith 
9774b9ad928SBarry Smith    Level: advanced
9784b9ad928SBarry Smith 
9794b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
9804b9ad928SBarry Smith 
98197177400SBarry Smith .seealso: PCMGSetLevels()
9824b9ad928SBarry Smith @*/
9837087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
9844b9ad928SBarry Smith {
985f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9864b9ad928SBarry Smith 
9874b9ad928SBarry Smith   PetscFunctionBegin;
9880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
989c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
99031567311SBarry Smith   mg->am = form;
9919dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
992c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
993c60c7ad4SBarry Smith   PetscFunctionReturn(0);
994c60c7ad4SBarry Smith }
995c60c7ad4SBarry Smith 
996c60c7ad4SBarry Smith /*@
997c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
998c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
999c60c7ad4SBarry Smith 
1000c60c7ad4SBarry Smith    Logically Collective on PC
1001c60c7ad4SBarry Smith 
1002c60c7ad4SBarry Smith    Input Parameter:
1003c60c7ad4SBarry Smith .  pc - the preconditioner context
1004c60c7ad4SBarry Smith 
1005c60c7ad4SBarry Smith    Output Parameter:
1006c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
1007c60c7ad4SBarry Smith 
1008c60c7ad4SBarry Smith 
1009c60c7ad4SBarry Smith    Level: advanced
1010c60c7ad4SBarry Smith 
1011c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
1012c60c7ad4SBarry Smith 
1013c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
1014c60c7ad4SBarry Smith @*/
1015c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1016c60c7ad4SBarry Smith {
1017c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1018c60c7ad4SBarry Smith 
1019c60c7ad4SBarry Smith   PetscFunctionBegin;
1020c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1021c60c7ad4SBarry Smith   *type = mg->am;
10224b9ad928SBarry Smith   PetscFunctionReturn(0);
10234b9ad928SBarry Smith }
10244b9ad928SBarry Smith 
10254b9ad928SBarry Smith /*@
10260d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
10274b9ad928SBarry Smith    complicated cycling.
10284b9ad928SBarry Smith 
1029ad4df100SBarry Smith    Logically Collective on PC
10304b9ad928SBarry Smith 
10314b9ad928SBarry Smith    Input Parameters:
1032c2be2410SBarry Smith +  pc - the multigrid context
1033c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
10344b9ad928SBarry Smith 
10354b9ad928SBarry Smith    Options Database Key:
1036c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
10374b9ad928SBarry Smith 
10384b9ad928SBarry Smith    Level: advanced
10394b9ad928SBarry Smith 
10404b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10414b9ad928SBarry Smith 
10420d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
10434b9ad928SBarry Smith @*/
10447087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
10454b9ad928SBarry Smith {
1046f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
1047f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
104879416396SBarry Smith   PetscInt     i,levels;
10494b9ad928SBarry Smith 
10504b9ad928SBarry Smith   PetscFunctionBegin;
10510700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
105269fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
10531a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1054f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10552fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
10564b9ad928SBarry Smith   PetscFunctionReturn(0);
10574b9ad928SBarry Smith }
10584b9ad928SBarry Smith 
10598cc2d5dfSBarry Smith /*@
10608cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
10618cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
10628cc2d5dfSBarry Smith 
1063ad4df100SBarry Smith    Logically Collective on PC
10648cc2d5dfSBarry Smith 
10658cc2d5dfSBarry Smith    Input Parameters:
10668cc2d5dfSBarry Smith +  pc - the multigrid context
10678cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
10688cc2d5dfSBarry Smith 
10698cc2d5dfSBarry Smith    Options Database Key:
1070e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
10718cc2d5dfSBarry Smith 
10728cc2d5dfSBarry Smith    Level: advanced
10738cc2d5dfSBarry Smith 
107495452b02SPatrick Sanan    Notes:
107595452b02SPatrick Sanan     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
10768cc2d5dfSBarry Smith 
10778cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10788cc2d5dfSBarry Smith 
10798cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
10808cc2d5dfSBarry Smith @*/
10817087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
10828cc2d5dfSBarry Smith {
1083f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
10848cc2d5dfSBarry Smith 
10858cc2d5dfSBarry Smith   PetscFunctionBegin;
10860700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1087c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
10882a21e185SBarry Smith   mg->cyclesperpcapply = n;
10898cc2d5dfSBarry Smith   PetscFunctionReturn(0);
10908cc2d5dfSBarry Smith }
10918cc2d5dfSBarry Smith 
10922134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1093fb15c04eSBarry Smith {
1094fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1095fb15c04eSBarry Smith 
1096fb15c04eSBarry Smith   PetscFunctionBegin;
10972134b1e4SBarry Smith   mg->galerkin = use;
1098fb15c04eSBarry Smith   PetscFunctionReturn(0);
1099fb15c04eSBarry Smith }
1100fb15c04eSBarry Smith 
1101c2be2410SBarry Smith /*@
110297177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
110382b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1104c2be2410SBarry Smith 
1105ad4df100SBarry Smith    Logically Collective on PC
1106c2be2410SBarry Smith 
1107c2be2410SBarry Smith    Input Parameters:
1108c91913d3SJed Brown +  pc - the multigrid context
11092134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1110c2be2410SBarry Smith 
1111c2be2410SBarry Smith    Options Database Key:
11122134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1113c2be2410SBarry Smith 
1114c2be2410SBarry Smith    Level: intermediate
1115c2be2410SBarry Smith 
111695452b02SPatrick Sanan    Notes:
111795452b02SPatrick Sanan     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
11181c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
11191c1aac46SBarry Smith 
1120c2be2410SBarry Smith .keywords: MG, set, Galerkin
1121c2be2410SBarry Smith 
11222134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
11233fc8bf9cSBarry Smith 
1124c2be2410SBarry Smith @*/
11252134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1126c2be2410SBarry Smith {
1127fb15c04eSBarry Smith   PetscErrorCode ierr;
1128c2be2410SBarry Smith 
1129c2be2410SBarry Smith   PetscFunctionBegin;
11300700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11312134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1132c2be2410SBarry Smith   PetscFunctionReturn(0);
1133c2be2410SBarry Smith }
1134c2be2410SBarry Smith 
11353fc8bf9cSBarry Smith /*@
11363fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
113782b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
11383fc8bf9cSBarry Smith 
11393fc8bf9cSBarry Smith    Not Collective
11403fc8bf9cSBarry Smith 
11413fc8bf9cSBarry Smith    Input Parameter:
11423fc8bf9cSBarry Smith .  pc - the multigrid context
11433fc8bf9cSBarry Smith 
11443fc8bf9cSBarry Smith    Output Parameter:
11452134b1e4SBarry 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
11463fc8bf9cSBarry Smith 
11473fc8bf9cSBarry Smith    Level: intermediate
11483fc8bf9cSBarry Smith 
11493fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
11503fc8bf9cSBarry Smith 
11512134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
11523fc8bf9cSBarry Smith 
11533fc8bf9cSBarry Smith @*/
11542134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
11553fc8bf9cSBarry Smith {
1156f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
11573fc8bf9cSBarry Smith 
11583fc8bf9cSBarry Smith   PetscFunctionBegin;
11590700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11602134b1e4SBarry Smith   *galerkin = mg->galerkin;
11613fc8bf9cSBarry Smith   PetscFunctionReturn(0);
11623fc8bf9cSBarry Smith }
11633fc8bf9cSBarry Smith 
11644b9ad928SBarry Smith /*@
116506792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1166710315b6SLawrence Mitchell    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1167710315b6SLawrence Mitchell    pre- and post-smoothing steps.
116806792cafSBarry Smith 
116906792cafSBarry Smith    Logically Collective on PC
117006792cafSBarry Smith 
117106792cafSBarry Smith    Input Parameters:
117206792cafSBarry Smith +  mg - the multigrid context
117306792cafSBarry Smith -  n - the number of smoothing steps
117406792cafSBarry Smith 
117506792cafSBarry Smith    Options Database Key:
117606792cafSBarry Smith +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
117706792cafSBarry Smith 
117806792cafSBarry Smith    Level: advanced
117906792cafSBarry Smith 
118095452b02SPatrick Sanan    Notes:
118195452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
118206792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
118306792cafSBarry Smith 
118406792cafSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
118506792cafSBarry Smith 
1186710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp()
118706792cafSBarry Smith @*/
118806792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
118906792cafSBarry Smith {
119006792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
119106792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
119206792cafSBarry Smith   PetscErrorCode ierr;
119306792cafSBarry Smith   PetscInt       i,levels;
119406792cafSBarry Smith 
119506792cafSBarry Smith   PetscFunctionBegin;
119606792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
119706792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
11981a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
119906792cafSBarry Smith   levels = mglevels[0]->levels;
120006792cafSBarry Smith 
120106792cafSBarry Smith   for (i=1; i<levels; i++) {
120206792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
120306792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
120406792cafSBarry Smith     mg->default_smoothu = n;
120506792cafSBarry Smith     mg->default_smoothd = n;
120606792cafSBarry Smith   }
120706792cafSBarry Smith   PetscFunctionReturn(0);
120806792cafSBarry Smith }
120906792cafSBarry Smith 
1210f442ab6aSBarry Smith /*@
1211f442ab6aSBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels
1212710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1213f442ab6aSBarry Smith 
1214f442ab6aSBarry Smith    Logically Collective on PC
1215f442ab6aSBarry Smith 
1216f442ab6aSBarry Smith    Input Parameters:
1217f442ab6aSBarry Smith .  pc - the preconditioner context
1218f442ab6aSBarry Smith 
1219f442ab6aSBarry Smith    Options Database Key:
1220f442ab6aSBarry Smith .  -pc_mg_distinct_smoothup
1221f442ab6aSBarry Smith 
1222f442ab6aSBarry Smith    Level: advanced
1223f442ab6aSBarry Smith 
122495452b02SPatrick Sanan    Notes:
122595452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
1226f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1227f442ab6aSBarry Smith 
1228f442ab6aSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1229f442ab6aSBarry Smith 
1230710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth()
1231f442ab6aSBarry Smith @*/
1232f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1233f442ab6aSBarry Smith {
1234f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1235f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1236f442ab6aSBarry Smith   PetscErrorCode ierr;
1237f442ab6aSBarry Smith   PetscInt       i,levels;
1238f442ab6aSBarry Smith   KSP            subksp;
1239f442ab6aSBarry Smith 
1240f442ab6aSBarry Smith   PetscFunctionBegin;
1241f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1242f442ab6aSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1243f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1244f442ab6aSBarry Smith 
1245f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1246710315b6SLawrence Mitchell     const char *prefix = NULL;
1247f442ab6aSBarry Smith     /* make sure smoother up and down are different */
1248f442ab6aSBarry Smith     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1249710315b6SLawrence Mitchell     ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr);
1250710315b6SLawrence Mitchell     ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr);
1251f442ab6aSBarry Smith     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1252f442ab6aSBarry Smith   }
1253f442ab6aSBarry Smith   PetscFunctionReturn(0);
1254f442ab6aSBarry Smith }
1255f442ab6aSBarry Smith 
12564b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
12574b9ad928SBarry Smith 
12583b09bd56SBarry Smith /*MC
1259ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
12603b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
12613b09bd56SBarry Smith 
12623b09bd56SBarry Smith    Options Database Keys:
12633b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1264391689abSStefano Zampini .  -pc_mg_cycle_type <v,w> - provide the cycle desired
12658c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
12663b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1267710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
12682134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
12698cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
12708cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1271e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
12728cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
12738cf6037eSBarry Smith                         to the binary output file called binaryoutput
12743b09bd56SBarry Smith 
127595452b02SPatrick Sanan    Notes:
12760b3c753eSRichard Tran Mills     If one uses a Krylov method such GMRES or CG as the smoother then one must use KSPFGMRES, KSPGCR, or KSPRICHARDSON as the outer Krylov method
12773b09bd56SBarry Smith 
12788cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
12798cf6037eSBarry Smith 
128023067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
128123067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
128223067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
128323067569SBarry Smith        residual is computed at the end of each cycle.
128423067569SBarry Smith 
12853b09bd56SBarry Smith    Level: intermediate
12863b09bd56SBarry Smith 
12878f87f92bSBarry Smith    Concepts: multigrid/multilevel
12883b09bd56SBarry Smith 
12898cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1290710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1291710315b6SLawrence Mitchell            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
129297177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12930d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12943b09bd56SBarry Smith M*/
12953b09bd56SBarry Smith 
12968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
12974b9ad928SBarry Smith {
1298f3fbd535SBarry Smith   PC_MG          *mg;
1299f3fbd535SBarry Smith   PetscErrorCode ierr;
1300f3fbd535SBarry Smith 
13014b9ad928SBarry Smith   PetscFunctionBegin;
1302b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1303f3fbd535SBarry Smith   pc->data     = (void*)mg;
1304f3fbd535SBarry Smith   mg->nlevels  = -1;
130510eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
13062134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1307f3fbd535SBarry Smith 
130837a44384SMark Adams   pc->useAmat = PETSC_TRUE;
130937a44384SMark Adams 
13104b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
13114b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1312a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
13134b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
13144b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
13154b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1316fb15c04eSBarry Smith 
1317fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
131897d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
131997d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
13204b9ad928SBarry Smith   PetscFunctionReturn(0);
13214b9ad928SBarry Smith }
1322