xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision c27ee7a386f5f8b417e4f67f23929a131c226f5e)
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;
154b9ad928SBarry Smith 
164b9ad928SBarry Smith   PetscFunctionBegin;
1763e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
1831567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
19c0decd05SBarry Smith   ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr);
2063e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2131567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2263e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2331567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2463e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
254b9ad928SBarry Smith 
264b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
2731567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
284b9ad928SBarry Smith       PetscReal rnorm;
2931567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
304b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3170441072SBarry Smith         if (rnorm < mg->abstol) {
324d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
3357622a8eSBarry 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);
344b9ad928SBarry Smith         } else {
354d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
3657622a8eSBarry 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);
374b9ad928SBarry Smith         }
384b9ad928SBarry Smith         PetscFunctionReturn(0);
394b9ad928SBarry Smith       }
404b9ad928SBarry Smith     }
414b9ad928SBarry Smith 
4231567311SBarry Smith     mgc = *(mglevelsin - 1);
4363e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4431567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
4563e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
46efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
474b9ad928SBarry Smith     while (cycles--) {
4831567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
494b9ad928SBarry Smith     }
5063e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5131567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5263e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5363e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5431567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
55c0decd05SBarry Smith     ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr);
5663e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith   }
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
61ace3abfcSBarry 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)
624b9ad928SBarry Smith {
63f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
64f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
65dfbe8321SBarry Smith   PetscErrorCode ierr;
66391689abSStefano Zampini   PC             tpc;
67391689abSStefano Zampini   PetscBool      changeu,changed;
68f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
694b9ad928SBarry Smith 
704b9ad928SBarry Smith   PetscFunctionBegin;
71a762d673SBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
72a762d673SBarry Smith   for (i=0; i<levels; i++) {
73a762d673SBarry Smith     if (!mglevels[i]->A) {
74a762d673SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
75a762d673SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
76a762d673SBarry Smith     }
77a762d673SBarry Smith   }
78391689abSStefano Zampini 
79391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
80391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
81391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
82391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
83391689abSStefano Zampini   if (!changed && !changeu) {
84391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
85f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
86391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
87391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
88391689abSStefano Zampini       Vec *vec;
89391689abSStefano Zampini 
90391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
91391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
927ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
93391689abSStefano Zampini     }
94391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
95391689abSStefano Zampini   }
96f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
974b9ad928SBarry Smith 
9831567311SBarry Smith   mg->rtol   = rtol;
9931567311SBarry Smith   mg->abstol = abstol;
10031567311SBarry Smith   mg->dtol   = dtol;
1014b9ad928SBarry Smith   if (rtol) {
1024b9ad928SBarry Smith     /* compute initial residual norm for relative convergence test */
1034b9ad928SBarry Smith     PetscReal rnorm;
1047319c654SBarry Smith     if (zeroguess) {
1057319c654SBarry Smith       ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr);
1067319c654SBarry Smith     } else {
107f3fbd535SBarry Smith       ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr);
1084b9ad928SBarry Smith       ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr);
1097319c654SBarry Smith     }
11031567311SBarry Smith     mg->ttol = PetscMax(rtol*rnorm,abstol);
1112fa5cd67SKarl Rupp   } else if (abstol) mg->ttol = abstol;
1122fa5cd67SKarl Rupp   else mg->ttol = 0.0;
1134b9ad928SBarry Smith 
1144d0a8057SBarry Smith   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
115336babb1SJed Brown      stop prematurely due to small residual */
1164d0a8057SBarry Smith   for (i=1; i<levels; i++) {
117f3fbd535SBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
118f3fbd535SBarry Smith     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
11923067569SBarry Smith       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
12023067569SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
121f3fbd535SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);
1224b9ad928SBarry Smith     }
1234d0a8057SBarry Smith   }
1244d0a8057SBarry Smith 
1254d0a8057SBarry Smith   *reason = (PCRichardsonConvergedReason)0;
1264d0a8057SBarry Smith   for (i=0; i<its; i++) {
127f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr);
1284d0a8057SBarry Smith     if (*reason) break;
1294d0a8057SBarry Smith   }
1304d0a8057SBarry Smith   if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
1314d0a8057SBarry Smith   *outits = i;
132391689abSStefano Zampini   if (!changed && !changeu) mglevels[levels-1]->b = NULL;
1334b9ad928SBarry Smith   PetscFunctionReturn(0);
1344b9ad928SBarry Smith }
1354b9ad928SBarry Smith 
1363aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc)
1373aeaf226SBarry Smith {
1383aeaf226SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1393aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1403aeaf226SBarry Smith   PetscErrorCode ierr;
1413aeaf226SBarry Smith   PetscInt       i,n;
1423aeaf226SBarry Smith 
1433aeaf226SBarry Smith   PetscFunctionBegin;
1443aeaf226SBarry Smith   if (mglevels) {
1453aeaf226SBarry Smith     n = mglevels[0]->levels;
1463aeaf226SBarry Smith     for (i=0; i<n-1; i++) {
1473aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr);
1483aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr);
1493aeaf226SBarry Smith       ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr);
1503aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr);
1513aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr);
1520de8493bSLawrence Mitchell       ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr);
15373250ac0SBarry Smith       ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr);
1543aeaf226SBarry Smith     }
155391689abSStefano Zampini     /* this is not null only if the smoother on the finest level
156391689abSStefano Zampini        changes the rhs during PreSolve */
157391689abSStefano Zampini     ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr);
1583aeaf226SBarry Smith 
1593aeaf226SBarry Smith     for (i=0; i<n; i++) {
1603aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1613aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1623aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1633aeaf226SBarry Smith       }
1643aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1653aeaf226SBarry Smith     }
1663aeaf226SBarry Smith   }
1673aeaf226SBarry Smith   PetscFunctionReturn(0);
1683aeaf226SBarry Smith }
1693aeaf226SBarry Smith 
17097d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
1714b9ad928SBarry Smith {
172dfbe8321SBarry Smith   PetscErrorCode ierr;
173f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
174ce94432eSBarry Smith   MPI_Comm       comm;
1753aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
17610eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
17710167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
178f3fbd535SBarry Smith   PetscInt       i;
179f3fbd535SBarry Smith   PetscMPIInt    size;
180f3fbd535SBarry Smith   const char     *prefix;
181f3fbd535SBarry Smith   PC             ipc;
1823aeaf226SBarry Smith   PetscInt       n;
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith   PetscFunctionBegin;
1850700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
186548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
187548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1881a97d7b6SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1893aeaf226SBarry Smith   if (mglevels) {
19010eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
1913aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
1923aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
1933aeaf226SBarry Smith     n    = mglevels[0]->levels;
1943aeaf226SBarry Smith     for (i=0; i<n; i++) {
1953aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1963aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
1973aeaf226SBarry Smith       }
1983aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
1993aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
2003aeaf226SBarry Smith     }
2013aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
2023aeaf226SBarry Smith   }
203f3fbd535SBarry Smith 
204f3fbd535SBarry Smith   mg->nlevels = levels;
205f3fbd535SBarry Smith 
206785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
2073bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
208f3fbd535SBarry Smith 
209f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
210f3fbd535SBarry Smith 
211a9db87a2SMatthew G Knepley   mg->stageApply = 0;
212f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
213b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
2142fa5cd67SKarl Rupp 
215f3fbd535SBarry Smith     mglevels[i]->level               = i;
216f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
21710eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
218336babb1SJed Brown     mg->default_smoothu              = 2;
219336babb1SJed Brown     mg->default_smoothd              = 2;
22063e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
22163e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
22263e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
22363e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
224f3fbd535SBarry Smith 
225f3fbd535SBarry Smith     if (comms) comm = comms[i];
226f3fbd535SBarry Smith     ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
227422a814eSBarry Smith     ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
2280ee9a668SBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
2290ee9a668SBarry Smith     ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
2300ee9a668SBarry Smith     ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
2310ee9a668SBarry Smith     if (i || levels == 1) {
2320ee9a668SBarry Smith       char tprefix[128];
2330ee9a668SBarry Smith 
234336babb1SJed Brown       ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2350059c7bdSJed Brown       ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
236336babb1SJed Brown       ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
237336babb1SJed Brown       ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
238336babb1SJed Brown       ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
2390ee9a668SBarry Smith       ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
240f3fbd535SBarry Smith 
2410ee9a668SBarry Smith       sprintf(tprefix,"mg_levels_%d_",(int)i);
2420ee9a668SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
2430ee9a668SBarry Smith     } else {
244f3fbd535SBarry Smith       ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
245f3fbd535SBarry Smith 
2469dbfc187SHong Zhang       /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
247f3fbd535SBarry Smith       ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
248f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
249f3fbd535SBarry Smith       ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
250f3fbd535SBarry Smith       if (size > 1) {
251f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
252f3fbd535SBarry Smith       } else {
253f3fbd535SBarry Smith         ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
254f3fbd535SBarry Smith       }
255753b7fb9SBarry Smith       ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
256f3fbd535SBarry Smith     }
2573bb1ff40SBarry Smith     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
2582fa5cd67SKarl Rupp 
259f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
26031567311SBarry Smith     mg->rtol             = 0.0;
26131567311SBarry Smith     mg->abstol           = 0.0;
26231567311SBarry Smith     mg->dtol             = 0.0;
26331567311SBarry Smith     mg->ttol             = 0.0;
26431567311SBarry Smith     mg->cyclesperpcapply = 1;
265f3fbd535SBarry Smith   }
266f3fbd535SBarry Smith   mg->levels = mglevels;
26710eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2684b9ad928SBarry Smith   PetscFunctionReturn(0);
2694b9ad928SBarry Smith }
2704b9ad928SBarry Smith 
27197d33e41SMatthew G. Knepley /*@C
27297d33e41SMatthew G. Knepley    PCMGSetLevels - Sets the number of levels to use with MG.
27397d33e41SMatthew G. Knepley    Must be called before any other MG routine.
27497d33e41SMatthew G. Knepley 
27597d33e41SMatthew G. Knepley    Logically Collective on PC
27697d33e41SMatthew G. Knepley 
27797d33e41SMatthew G. Knepley    Input Parameters:
27897d33e41SMatthew G. Knepley +  pc - the preconditioner context
27997d33e41SMatthew G. Knepley .  levels - the number of levels
28097d33e41SMatthew G. Knepley -  comms - optional communicators for each level; this is to allow solving the coarser problems
28197d33e41SMatthew G. Knepley            on smaller sets of processors.
28297d33e41SMatthew G. Knepley 
28397d33e41SMatthew G. Knepley    Level: intermediate
28497d33e41SMatthew G. Knepley 
28597d33e41SMatthew G. Knepley    Notes:
28697d33e41SMatthew G. Knepley      If the number of levels is one then the multigrid uses the -mg_levels prefix
28797d33e41SMatthew G. Knepley   for setting the level options rather than the -mg_coarse prefix.
28897d33e41SMatthew G. Knepley 
28997d33e41SMatthew G. Knepley .keywords: MG, set, levels, multigrid
29097d33e41SMatthew G. Knepley 
29197d33e41SMatthew G. Knepley .seealso: PCMGSetType(), PCMGGetLevels()
29297d33e41SMatthew G. Knepley @*/
29397d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
29497d33e41SMatthew G. Knepley {
29597d33e41SMatthew G. Knepley   PetscErrorCode ierr;
29697d33e41SMatthew G. Knepley 
29797d33e41SMatthew G. Knepley   PetscFunctionBegin;
29897d33e41SMatthew G. Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
29997d33e41SMatthew G. Knepley   if (comms) PetscValidPointer(comms,3);
30037b5128cSMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr);
30197d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
30297d33e41SMatthew G. Knepley }
30397d33e41SMatthew G. Knepley 
304c07bf074SBarry Smith 
305c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
306c07bf074SBarry Smith {
307c07bf074SBarry Smith   PetscErrorCode ierr;
308a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
309a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
310a06653b4SBarry Smith   PetscInt       i,n;
311c07bf074SBarry Smith 
312c07bf074SBarry Smith   PetscFunctionBegin;
313a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
314a06653b4SBarry Smith   if (mglevels) {
315a06653b4SBarry Smith     n = mglevels[0]->levels;
316a06653b4SBarry Smith     for (i=0; i<n; i++) {
317a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
3186bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
319a06653b4SBarry Smith       }
3206bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
321a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
322a06653b4SBarry Smith     }
323a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
324a06653b4SBarry Smith   }
325c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
326f3fbd535SBarry Smith   PetscFunctionReturn(0);
327f3fbd535SBarry Smith }
328f3fbd535SBarry Smith 
329f3fbd535SBarry Smith 
330f3fbd535SBarry Smith 
33109573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
33209573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
33309573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
334f3fbd535SBarry Smith 
335f3fbd535SBarry Smith /*
336f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
337f3fbd535SBarry Smith              or full cycle of multigrid.
338f3fbd535SBarry Smith 
339f3fbd535SBarry Smith   Note:
340f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
341f3fbd535SBarry Smith */
342f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
343f3fbd535SBarry Smith {
344f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
345f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
346f3fbd535SBarry Smith   PetscErrorCode ierr;
347391689abSStefano Zampini   PC             tpc;
348f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
349391689abSStefano Zampini   PetscBool      changeu,changed;
350f3fbd535SBarry Smith 
351f3fbd535SBarry Smith   PetscFunctionBegin;
352a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
353e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
354298cc208SBarry Smith   for (i=0; i<levels; i++) {
355e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
35623ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
357298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
358e1d8e5deSBarry Smith     }
359e1d8e5deSBarry Smith   }
360e1d8e5deSBarry Smith 
361391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
362391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
363391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
364391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
365391689abSStefano Zampini   if (!changeu && !changed) {
366391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
367f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
368391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
369391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
370391689abSStefano Zampini       Vec *vec;
371391689abSStefano Zampini 
372391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
373391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
3747ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
375391689abSStefano Zampini     }
376391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
377391689abSStefano Zampini   }
378f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
379391689abSStefano Zampini 
38031567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
381f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
38231567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
3830298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
384f3fbd535SBarry Smith     }
3852fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
38631567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
3872fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
38831567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
3892fa5cd67SKarl Rupp   } else {
390f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
391f3fbd535SBarry Smith   }
392a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
393391689abSStefano Zampini   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
394f3fbd535SBarry Smith   PetscFunctionReturn(0);
395f3fbd535SBarry Smith }
396f3fbd535SBarry Smith 
397f3fbd535SBarry Smith 
3984416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
399f3fbd535SBarry Smith {
400f3fbd535SBarry Smith   PetscErrorCode   ierr;
401710315b6SLawrence Mitchell   PetscInt         levels,cycles;
4022134b1e4SBarry Smith   PetscBool        flg;
403f3fbd535SBarry Smith   PC_MG            *mg = (PC_MG*)pc->data;
4043d3eaba7SBarry Smith   PC_MG_Levels     **mglevels;
405f3fbd535SBarry Smith   PCMGType         mgtype;
406f3fbd535SBarry Smith   PCMGCycleType    mgctype;
4072134b1e4SBarry Smith   PCMGGalerkinType gtype;
408f3fbd535SBarry Smith 
409f3fbd535SBarry Smith   PetscFunctionBegin;
4101d743356SStefano Zampini   levels = PetscMax(mg->nlevels,1);
411e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
412f3fbd535SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
4131a97d7b6SStefano Zampini   if (!flg && !mg->levels && pc->dm) {
414eb3f98d2SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
415eb3f98d2SBarry Smith     levels++;
4163aeaf226SBarry Smith     mg->usedmfornumberoflevels = PETSC_TRUE;
417eb3f98d2SBarry Smith   }
4180298fd71SBarry Smith   ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
4193aeaf226SBarry Smith   mglevels = mg->levels;
4203aeaf226SBarry Smith 
421f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
422f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
423f3fbd535SBarry Smith   if (flg) {
424f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
4252fa5cd67SKarl Rupp   }
4262134b1e4SBarry Smith   gtype = mg->galerkin;
4272134b1e4SBarry Smith   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
4282134b1e4SBarry Smith   if (flg) {
4292134b1e4SBarry Smith     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
430f3fbd535SBarry Smith   }
431f56b1265SBarry Smith   flg = PETSC_FALSE;
432f442ab6aSBarry Smith   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
433f442ab6aSBarry Smith   if (flg) {
434f442ab6aSBarry Smith     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
435f442ab6aSBarry Smith   }
43631567311SBarry Smith   mgtype = mg->am;
437f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
438f3fbd535SBarry Smith   if (flg) {
439f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
440f3fbd535SBarry Smith   }
44131567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
4425363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
443f3fbd535SBarry Smith     if (flg) {
444f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
445f3fbd535SBarry Smith     }
446f3fbd535SBarry Smith   }
447f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4480298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
449f3fbd535SBarry Smith   if (flg) {
450f3fbd535SBarry Smith     PetscInt i;
451f3fbd535SBarry Smith     char     eventname[128];
4521a97d7b6SStefano Zampini 
453f3fbd535SBarry Smith     levels = mglevels[0]->levels;
454f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
455f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
45663e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
457f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
45863e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
459f3fbd535SBarry Smith       if (i) {
460f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
46163e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
462f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
46363e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
464f3fbd535SBarry Smith       }
465f3fbd535SBarry Smith     }
466eec35531SMatthew G Knepley 
467ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
468eec35531SMatthew G Knepley     {
469eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
470eec35531SMatthew G Knepley       PetscStageLog stageLog;
471eec35531SMatthew G Knepley       PetscInt      st;
472eec35531SMatthew G Knepley 
473eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
474eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
475eec35531SMatthew G Knepley         PetscBool same;
476eec35531SMatthew G Knepley 
477eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
4782fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
479eec35531SMatthew G Knepley       }
480eec35531SMatthew G Knepley       if (!mg->stageApply) {
481eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
482eec35531SMatthew G Knepley       }
483eec35531SMatthew G Knepley     }
484ec60ca73SSatish Balay #endif
485f3fbd535SBarry Smith   }
486f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
487f3fbd535SBarry Smith   PetscFunctionReturn(0);
488f3fbd535SBarry Smith }
489f3fbd535SBarry Smith 
4906a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
4916a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
4922134b1e4SBarry Smith const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0};
493f3fbd535SBarry Smith 
4949804daf3SBarry Smith #include <petscdraw.h>
495f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
496f3fbd535SBarry Smith {
497f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
498f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
499f3fbd535SBarry Smith   PetscErrorCode ierr;
500e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
501219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
502f3fbd535SBarry Smith 
503f3fbd535SBarry Smith   PetscFunctionBegin;
504251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5055b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
506219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
507f3fbd535SBarry Smith   if (iascii) {
508e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
509efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
51031567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
51131567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
512f3fbd535SBarry Smith     }
5132134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
514f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
5152134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
5162134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
5172134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
5182134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
5192134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
5202134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
5214f66f45eSBarry Smith     } else {
5224f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
523f3fbd535SBarry Smith     }
5245adeb434SBarry Smith     if (mg->view){
5255adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
5265adeb434SBarry Smith     }
527f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
528f3fbd535SBarry Smith       if (!i) {
52989cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
530f3fbd535SBarry Smith       } else {
53189cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
532f3fbd535SBarry Smith       }
533f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
534f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
535f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
536f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
537f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
538f3fbd535SBarry Smith       } else if (i) {
53989cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
540f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
541f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
542f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
543f3fbd535SBarry Smith       }
544f3fbd535SBarry Smith     }
5455b0b0462SBarry Smith   } else if (isbinary) {
5465b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5475b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5485b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5495b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5505b0b0462SBarry Smith       }
5515b0b0462SBarry Smith     }
552219b1045SBarry Smith   } else if (isdraw) {
553219b1045SBarry Smith     PetscDraw draw;
554b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
555219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
556219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5570298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
558219b1045SBarry Smith     bottom = y - th;
559219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
560b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
561219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
562219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
563219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
564b4375e8dSPeter Brune       } else {
565b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
566b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
567b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
568b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
569b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
570b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
571b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
572b4375e8dSPeter Brune       }
5730298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
5741cd381d6SBarry Smith       bottom -= th;
575219b1045SBarry Smith     }
5765b0b0462SBarry Smith   }
577f3fbd535SBarry Smith   PetscFunctionReturn(0);
578f3fbd535SBarry Smith }
579f3fbd535SBarry Smith 
580af0996ceSBarry Smith #include <petsc/private/dmimpl.h>
581af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
582cab2e9ccSBarry Smith 
583f3fbd535SBarry Smith /*
584f3fbd535SBarry Smith     Calls setup for the KSP on each level
585f3fbd535SBarry Smith */
586f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
587f3fbd535SBarry Smith {
588f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
589f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
590f3fbd535SBarry Smith   PetscErrorCode ierr;
5911a97d7b6SStefano Zampini   PetscInt       i,n;
59298e04984SBarry Smith   PC             cpc;
5933db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
594f3fbd535SBarry Smith   Mat            dA,dB;
595f3fbd535SBarry Smith   Vec            tvec;
596218a07d4SBarry Smith   DM             *dms;
597649052a6SBarry Smith   PetscViewer    viewer = 0;
5982134b1e4SBarry Smith   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
599f3fbd535SBarry Smith 
600f3fbd535SBarry Smith   PetscFunctionBegin;
6011a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
6021a97d7b6SStefano Zampini   n = mglevels[0]->levels;
60301bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
6043aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
6053aeaf226SBarry Smith     PetscInt levels;
6063aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
6073aeaf226SBarry Smith     levels++;
6083aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
6090298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
6103aeaf226SBarry Smith       n        = levels;
6113aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
6123aeaf226SBarry Smith       mglevels = mg->levels;
6133aeaf226SBarry Smith     }
6143aeaf226SBarry Smith   }
61598e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
6163aeaf226SBarry Smith 
617f3fbd535SBarry Smith 
618f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
619f3fbd535SBarry Smith   /* so use those from global PC */
620f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
6210298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
62298e04984SBarry Smith   if (opsset) {
62398e04984SBarry Smith     Mat mmat;
62423ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
62598e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
62698e04984SBarry Smith   }
6274a5f07a5SMark F. Adams 
6284a5f07a5SMark F. Adams   if (!opsset) {
62971b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
6304a5f07a5SMark F. Adams     if (use_amat) {
631f3fbd535SBarry 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);
63223ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
63369858f1bSStefano Zampini     } else {
6344a5f07a5SMark 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);
63523ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
6364a5f07a5SMark F. Adams     }
6374a5f07a5SMark F. Adams   }
638f3fbd535SBarry Smith 
6399c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
6409c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
6419c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
6429c7ffb39SBarry Smith       continue;
6439c7ffb39SBarry Smith     }
6449c7ffb39SBarry Smith   }
6452134b1e4SBarry Smith 
6462134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
6472134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
6482134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
6492134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
6502134b1e4SBarry Smith   }
6512134b1e4SBarry Smith 
6522134b1e4SBarry Smith 
6539c7ffb39SBarry Smith   /*
6549c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6559c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6569c7ffb39SBarry Smith   */
6572134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
6582d2b81a6SBarry Smith 	/* construct the interpolation from the DMs */
6592e499ae9SBarry Smith     Mat p;
66073250ac0SBarry Smith     Vec rscale;
661785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
662218a07d4SBarry Smith     dms[n-1] = pc->dm;
663ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
664ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
66541fce8fdSHong Zhang 	/*
66641fce8fdSHong Zhang 	   Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
66741fce8fdSHong Zhang 	   Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
66841fce8fdSHong Zhang 	   But it is safe to use -dm_mat_type.
66941fce8fdSHong Zhang 
67041fce8fdSHong Zhang 	   The mat type should not be hardcoded like this, we need to find a better way.
67141fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
67241fce8fdSHong Zhang     */
673218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
674942e3340SBarry Smith       DMKSP     kdm;
675eab52d2bSLawrence Mitchell       PetscBool dmhasrestrict, dmhasinject;
6763c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
6772134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
678*c27ee7a3SPatrick Farrell       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
679*c27ee7a3SPatrick Farrell         ierr = KSPSetDM(mglevels[i]->smoothu,dms[i]);CHKERRQ(ierr);
680*c27ee7a3SPatrick Farrell         if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);CHKERRQ(ierr);}
681*c27ee7a3SPatrick Farrell       }
682942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
683d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
684d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
6850298fd71SBarry Smith       kdm->ops->computerhs = NULL;
6860298fd71SBarry Smith       kdm->rhsctx          = NULL;
68724c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
688e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
6892d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
69000ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
69173250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
6926bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
693218a07d4SBarry Smith       }
6943ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
6953ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
6963ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
6973ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
6983ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
6993ad4599aSBarry Smith       }
700eab52d2bSLawrence Mitchell       ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr);
701eab52d2bSLawrence Mitchell       if (dmhasinject && !mglevels[i+1]->inject){
702eab52d2bSLawrence Mitchell         ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr);
703eab52d2bSLawrence Mitchell         ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr);
704eab52d2bSLawrence Mitchell         ierr = MatDestroy(&p);CHKERRQ(ierr);
705eab52d2bSLawrence Mitchell       }
70624c3aa18SBarry Smith     }
7072d2b81a6SBarry Smith 
708ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
7092d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
710ce4cda84SJed Brown   }
711cab2e9ccSBarry Smith 
712ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
7132134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
714cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
715cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
716*c27ee7a3SPatrick Farrell     if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) {
717*c27ee7a3SPatrick Farrell       ierr = KSPSetDM(mglevels[n-1]->smoothu,pc->dm);CHKERRQ(ierr);
718*c27ee7a3SPatrick Farrell       ierr = KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);CHKERRQ(ierr);
719*c27ee7a3SPatrick Farrell     }
720218a07d4SBarry Smith   }
721218a07d4SBarry Smith 
7222134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
7232134b1e4SBarry Smith     Mat       A,B;
7242134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
7252134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
7262134b1e4SBarry Smith 
7272134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
7282134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
7292134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
730f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
731b9367914SBarry 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");
732b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
733b9367914SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
734b9367914SBarry Smith       }
735b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
736b9367914SBarry Smith         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
737b9367914SBarry Smith       }
7382134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
7392134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
7402134b1e4SBarry Smith       }
7412134b1e4SBarry Smith       if (doA) {
7422df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
7432134b1e4SBarry Smith       }
7442134b1e4SBarry Smith       if (doB) {
7452df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
746b9367914SBarry Smith       }
7472134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
7482134b1e4SBarry Smith       if (!doA && dAeqdB) {
7492134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
7502134b1e4SBarry Smith         A = B;
7512134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX ) {
7522134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
7532134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
754b9367914SBarry Smith       }
7552134b1e4SBarry Smith       if (!doB && dAeqdB) {
7562134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
7572134b1e4SBarry Smith         B = A;
7582134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
75923ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
7602134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
7617d537d92SJed Brown       }
7622134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
7632134b1e4SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
7642134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
7652134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
7662134b1e4SBarry Smith       }
7672134b1e4SBarry Smith       dA = A;
768f3fbd535SBarry Smith       dB = B;
769f3fbd535SBarry Smith     }
770f3fbd535SBarry Smith   }
7712134b1e4SBarry Smith   if (needRestricts && pc->dm && pc->dm->x) {
772cab2e9ccSBarry Smith     /* need to restrict Jacobian location to coarser meshes for evaluation */
773cab2e9ccSBarry Smith     for (i=n-2; i>-1; i--) {
774c88c5224SJed Brown       Mat R;
775c88c5224SJed Brown       Vec rscale;
776cab2e9ccSBarry Smith       if (!mglevels[i]->smoothd->dm->x) {
777cab2e9ccSBarry Smith         Vec *vecs;
7782a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr);
779cab2e9ccSBarry Smith         mglevels[i]->smoothd->dm->x = vecs[0];
780cab2e9ccSBarry Smith         ierr = PetscFree(vecs);CHKERRQ(ierr);
781cab2e9ccSBarry Smith       }
782c88c5224SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr);
783c88c5224SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
784c88c5224SJed Brown       ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr);
785c88c5224SJed Brown       ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr);
786cab2e9ccSBarry Smith     }
787f3fbd535SBarry Smith   }
7882134b1e4SBarry Smith   if (needRestricts && pc->dm) {
789caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
790ccceb50cSJed Brown       DM  dmfine,dmcoarse;
791caa4e7f2SJed Brown       Mat Restrict,Inject;
792caa4e7f2SJed Brown       Vec rscale;
793ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
794ccceb50cSJed Brown       ierr   = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
795caa4e7f2SJed Brown       ierr   = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
796caa4e7f2SJed Brown       ierr   = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
797eab52d2bSLawrence Mitchell       ierr   = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr);
798caa4e7f2SJed Brown       ierr   = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
799caa4e7f2SJed Brown     }
800caa4e7f2SJed Brown   }
801f3fbd535SBarry Smith 
802f3fbd535SBarry Smith   if (!pc->setupcalled) {
803f3fbd535SBarry Smith     for (i=0; i<n; i++) {
804f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
805f3fbd535SBarry Smith     }
806f3fbd535SBarry Smith     for (i=1; i<n; i++) {
807f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
808f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
809f3fbd535SBarry Smith       }
810f3fbd535SBarry Smith     }
8113ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
812f3fbd535SBarry Smith     for (i=1; i<n; i++) {
8133ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
8143ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
815f3fbd535SBarry Smith     }
816f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
817f3fbd535SBarry Smith       if (!mglevels[i]->b) {
818f3fbd535SBarry Smith         Vec *vec;
8192a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
820f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
8216bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
822f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
823f3fbd535SBarry Smith       }
824f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
825f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
826f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
8276bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
828f3fbd535SBarry Smith       }
829f3fbd535SBarry Smith       if (!mglevels[i]->x) {
830f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
831f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
8326bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
833f3fbd535SBarry Smith       }
834f3fbd535SBarry Smith     }
835f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
836f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
837f3fbd535SBarry Smith       Vec *vec;
8382a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
839f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
8406bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
841f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
842f3fbd535SBarry Smith     }
843f3fbd535SBarry Smith   }
844f3fbd535SBarry Smith 
84598e04984SBarry Smith   if (pc->dm) {
84698e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
84798e04984SBarry Smith     for (i=0; i<n-1; i++) {
84898e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
84998e04984SBarry Smith     }
85098e04984SBarry Smith   }
851f3fbd535SBarry Smith 
852f3fbd535SBarry Smith   for (i=1; i<n; i++) {
8532a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
854f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
855f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
856f3fbd535SBarry Smith     }
85763e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
858f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
859c0decd05SBarry Smith     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
860899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
861899639b0SHong Zhang     }
86263e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
863d42688cbSBarry Smith     if (!mglevels[i]->residual) {
864d42688cbSBarry Smith       Mat mat;
86513842ffbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
86654b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
867d42688cbSBarry Smith     }
868f3fbd535SBarry Smith   }
869f3fbd535SBarry Smith   for (i=1; i<n; i++) {
870f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
871f3fbd535SBarry Smith       Mat downmat,downpmat;
872f3fbd535SBarry Smith 
873f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
8740298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
875f3fbd535SBarry Smith       if (!opsset) {
87623ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
87723ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
878f3fbd535SBarry Smith       }
879f3fbd535SBarry Smith 
880f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
88163e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
882f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
883c0decd05SBarry Smith       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
884899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
885899639b0SHong Zhang       }
88663e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
887f3fbd535SBarry Smith     }
888f3fbd535SBarry Smith   }
889f3fbd535SBarry Smith 
89063e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
891f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
892c0decd05SBarry Smith   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
893899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
894899639b0SHong Zhang   }
89563e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
896f3fbd535SBarry Smith 
897f3fbd535SBarry Smith   /*
898f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
899e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
900f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
901f3fbd535SBarry Smith 
902f3fbd535SBarry Smith    Only support one or the other at the same time.
903f3fbd535SBarry Smith   */
904f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
905c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
906ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
907f3fbd535SBarry Smith   dump = PETSC_FALSE;
908f3fbd535SBarry Smith #endif
909c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
910ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
911f3fbd535SBarry Smith 
912f3fbd535SBarry Smith   if (viewer) {
913f3fbd535SBarry Smith     for (i=1; i<n; i++) {
914f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
915f3fbd535SBarry Smith     }
916f3fbd535SBarry Smith     for (i=0; i<n; i++) {
917f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
918f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
919f3fbd535SBarry Smith     }
920f3fbd535SBarry Smith   }
921f3fbd535SBarry Smith   PetscFunctionReturn(0);
922f3fbd535SBarry Smith }
923f3fbd535SBarry Smith 
924f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
925f3fbd535SBarry Smith 
92697d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
92797d33e41SMatthew G. Knepley {
92897d33e41SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
92997d33e41SMatthew G. Knepley 
93097d33e41SMatthew G. Knepley   PetscFunctionBegin;
93197d33e41SMatthew G. Knepley   *levels = mg->nlevels;
93297d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
93397d33e41SMatthew G. Knepley }
93497d33e41SMatthew G. Knepley 
9354b9ad928SBarry Smith /*@
93697177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
9374b9ad928SBarry Smith 
9384b9ad928SBarry Smith    Not Collective
9394b9ad928SBarry Smith 
9404b9ad928SBarry Smith    Input Parameter:
9414b9ad928SBarry Smith .  pc - the preconditioner context
9424b9ad928SBarry Smith 
9434b9ad928SBarry Smith    Output parameter:
9444b9ad928SBarry Smith .  levels - the number of levels
9454b9ad928SBarry Smith 
9464b9ad928SBarry Smith    Level: advanced
9474b9ad928SBarry Smith 
9484b9ad928SBarry Smith .keywords: MG, get, levels, multigrid
9494b9ad928SBarry Smith 
95097177400SBarry Smith .seealso: PCMGSetLevels()
9514b9ad928SBarry Smith @*/
9527087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
9534b9ad928SBarry Smith {
954e952c529SMatthew G. Knepley   PetscErrorCode ierr;
9554b9ad928SBarry Smith 
9564b9ad928SBarry Smith   PetscFunctionBegin;
9570700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9584482741eSBarry Smith   PetscValidIntPointer(levels,2);
95997d33e41SMatthew G. Knepley   *levels = 0;
96037b5128cSMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr);
9614b9ad928SBarry Smith   PetscFunctionReturn(0);
9624b9ad928SBarry Smith }
9634b9ad928SBarry Smith 
9644b9ad928SBarry Smith /*@
96597177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
9664b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
9674b9ad928SBarry Smith 
968ad4df100SBarry Smith    Logically Collective on PC
9694b9ad928SBarry Smith 
9704b9ad928SBarry Smith    Input Parameters:
9714b9ad928SBarry Smith +  pc - the preconditioner context
9729dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
9739dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
9744b9ad928SBarry Smith 
9754b9ad928SBarry Smith    Options Database Key:
9764b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
9774b9ad928SBarry Smith    additive, full, kaskade
9784b9ad928SBarry Smith 
9794b9ad928SBarry Smith    Level: advanced
9804b9ad928SBarry Smith 
9814b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
9824b9ad928SBarry Smith 
98397177400SBarry Smith .seealso: PCMGSetLevels()
9844b9ad928SBarry Smith @*/
9857087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
9864b9ad928SBarry Smith {
987f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
9884b9ad928SBarry Smith 
9894b9ad928SBarry Smith   PetscFunctionBegin;
9900700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
991c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
99231567311SBarry Smith   mg->am = form;
9939dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
994c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
995c60c7ad4SBarry Smith   PetscFunctionReturn(0);
996c60c7ad4SBarry Smith }
997c60c7ad4SBarry Smith 
998c60c7ad4SBarry Smith /*@
999c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
1000c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
1001c60c7ad4SBarry Smith 
1002c60c7ad4SBarry Smith    Logically Collective on PC
1003c60c7ad4SBarry Smith 
1004c60c7ad4SBarry Smith    Input Parameter:
1005c60c7ad4SBarry Smith .  pc - the preconditioner context
1006c60c7ad4SBarry Smith 
1007c60c7ad4SBarry Smith    Output Parameter:
1008c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
1009c60c7ad4SBarry Smith 
1010c60c7ad4SBarry Smith 
1011c60c7ad4SBarry Smith    Level: advanced
1012c60c7ad4SBarry Smith 
1013c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
1014c60c7ad4SBarry Smith 
1015c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
1016c60c7ad4SBarry Smith @*/
1017c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1018c60c7ad4SBarry Smith {
1019c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1020c60c7ad4SBarry Smith 
1021c60c7ad4SBarry Smith   PetscFunctionBegin;
1022c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1023c60c7ad4SBarry Smith   *type = mg->am;
10244b9ad928SBarry Smith   PetscFunctionReturn(0);
10254b9ad928SBarry Smith }
10264b9ad928SBarry Smith 
10274b9ad928SBarry Smith /*@
10280d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
10294b9ad928SBarry Smith    complicated cycling.
10304b9ad928SBarry Smith 
1031ad4df100SBarry Smith    Logically Collective on PC
10324b9ad928SBarry Smith 
10334b9ad928SBarry Smith    Input Parameters:
1034c2be2410SBarry Smith +  pc - the multigrid context
1035c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
10364b9ad928SBarry Smith 
10374b9ad928SBarry Smith    Options Database Key:
1038c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
10394b9ad928SBarry Smith 
10404b9ad928SBarry Smith    Level: advanced
10414b9ad928SBarry Smith 
10424b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10434b9ad928SBarry Smith 
10440d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
10454b9ad928SBarry Smith @*/
10467087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
10474b9ad928SBarry Smith {
1048f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
1049f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
105079416396SBarry Smith   PetscInt     i,levels;
10514b9ad928SBarry Smith 
10524b9ad928SBarry Smith   PetscFunctionBegin;
10530700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
105469fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
10551a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1056f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10572fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
10584b9ad928SBarry Smith   PetscFunctionReturn(0);
10594b9ad928SBarry Smith }
10604b9ad928SBarry Smith 
10618cc2d5dfSBarry Smith /*@
10628cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
10638cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
10648cc2d5dfSBarry Smith 
1065ad4df100SBarry Smith    Logically Collective on PC
10668cc2d5dfSBarry Smith 
10678cc2d5dfSBarry Smith    Input Parameters:
10688cc2d5dfSBarry Smith +  pc - the multigrid context
10698cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
10708cc2d5dfSBarry Smith 
10718cc2d5dfSBarry Smith    Options Database Key:
1072e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
10738cc2d5dfSBarry Smith 
10748cc2d5dfSBarry Smith    Level: advanced
10758cc2d5dfSBarry Smith 
107695452b02SPatrick Sanan    Notes:
107795452b02SPatrick Sanan     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
10788cc2d5dfSBarry Smith 
10798cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
10808cc2d5dfSBarry Smith 
10818cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
10828cc2d5dfSBarry Smith @*/
10837087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
10848cc2d5dfSBarry Smith {
1085f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
10868cc2d5dfSBarry Smith 
10878cc2d5dfSBarry Smith   PetscFunctionBegin;
10880700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1089c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
10902a21e185SBarry Smith   mg->cyclesperpcapply = n;
10918cc2d5dfSBarry Smith   PetscFunctionReturn(0);
10928cc2d5dfSBarry Smith }
10938cc2d5dfSBarry Smith 
10942134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1095fb15c04eSBarry Smith {
1096fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1097fb15c04eSBarry Smith 
1098fb15c04eSBarry Smith   PetscFunctionBegin;
10992134b1e4SBarry Smith   mg->galerkin = use;
1100fb15c04eSBarry Smith   PetscFunctionReturn(0);
1101fb15c04eSBarry Smith }
1102fb15c04eSBarry Smith 
1103c2be2410SBarry Smith /*@
110497177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
110582b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1106c2be2410SBarry Smith 
1107ad4df100SBarry Smith    Logically Collective on PC
1108c2be2410SBarry Smith 
1109c2be2410SBarry Smith    Input Parameters:
1110c91913d3SJed Brown +  pc - the multigrid context
11112134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1112c2be2410SBarry Smith 
1113c2be2410SBarry Smith    Options Database Key:
11142134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1115c2be2410SBarry Smith 
1116c2be2410SBarry Smith    Level: intermediate
1117c2be2410SBarry Smith 
111895452b02SPatrick Sanan    Notes:
111995452b02SPatrick Sanan     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
11201c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
11211c1aac46SBarry Smith 
1122c2be2410SBarry Smith .keywords: MG, set, Galerkin
1123c2be2410SBarry Smith 
11242134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
11253fc8bf9cSBarry Smith 
1126c2be2410SBarry Smith @*/
11272134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1128c2be2410SBarry Smith {
1129fb15c04eSBarry Smith   PetscErrorCode ierr;
1130c2be2410SBarry Smith 
1131c2be2410SBarry Smith   PetscFunctionBegin;
11320700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11332134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1134c2be2410SBarry Smith   PetscFunctionReturn(0);
1135c2be2410SBarry Smith }
1136c2be2410SBarry Smith 
11373fc8bf9cSBarry Smith /*@
11383fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
113982b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
11403fc8bf9cSBarry Smith 
11413fc8bf9cSBarry Smith    Not Collective
11423fc8bf9cSBarry Smith 
11433fc8bf9cSBarry Smith    Input Parameter:
11443fc8bf9cSBarry Smith .  pc - the multigrid context
11453fc8bf9cSBarry Smith 
11463fc8bf9cSBarry Smith    Output Parameter:
11472134b1e4SBarry 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
11483fc8bf9cSBarry Smith 
11493fc8bf9cSBarry Smith    Level: intermediate
11503fc8bf9cSBarry Smith 
11513fc8bf9cSBarry Smith .keywords: MG, set, Galerkin
11523fc8bf9cSBarry Smith 
11532134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
11543fc8bf9cSBarry Smith 
11553fc8bf9cSBarry Smith @*/
11562134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
11573fc8bf9cSBarry Smith {
1158f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
11593fc8bf9cSBarry Smith 
11603fc8bf9cSBarry Smith   PetscFunctionBegin;
11610700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11622134b1e4SBarry Smith   *galerkin = mg->galerkin;
11633fc8bf9cSBarry Smith   PetscFunctionReturn(0);
11643fc8bf9cSBarry Smith }
11653fc8bf9cSBarry Smith 
11664b9ad928SBarry Smith /*@
116706792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1168710315b6SLawrence Mitchell    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1169710315b6SLawrence Mitchell    pre- and post-smoothing steps.
117006792cafSBarry Smith 
117106792cafSBarry Smith    Logically Collective on PC
117206792cafSBarry Smith 
117306792cafSBarry Smith    Input Parameters:
117406792cafSBarry Smith +  mg - the multigrid context
117506792cafSBarry Smith -  n - the number of smoothing steps
117606792cafSBarry Smith 
117706792cafSBarry Smith    Options Database Key:
117806792cafSBarry Smith +  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
117906792cafSBarry Smith 
118006792cafSBarry Smith    Level: advanced
118106792cafSBarry Smith 
118295452b02SPatrick Sanan    Notes:
118395452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
118406792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
118506792cafSBarry Smith 
118606792cafSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
118706792cafSBarry Smith 
1188710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp()
118906792cafSBarry Smith @*/
119006792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
119106792cafSBarry Smith {
119206792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
119306792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
119406792cafSBarry Smith   PetscErrorCode ierr;
119506792cafSBarry Smith   PetscInt       i,levels;
119606792cafSBarry Smith 
119706792cafSBarry Smith   PetscFunctionBegin;
119806792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
119906792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
12001a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
120106792cafSBarry Smith   levels = mglevels[0]->levels;
120206792cafSBarry Smith 
120306792cafSBarry Smith   for (i=1; i<levels; i++) {
120406792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
120506792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
120606792cafSBarry Smith     mg->default_smoothu = n;
120706792cafSBarry Smith     mg->default_smoothd = n;
120806792cafSBarry Smith   }
120906792cafSBarry Smith   PetscFunctionReturn(0);
121006792cafSBarry Smith }
121106792cafSBarry Smith 
1212f442ab6aSBarry Smith /*@
1213f442ab6aSBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels
1214710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1215f442ab6aSBarry Smith 
1216f442ab6aSBarry Smith    Logically Collective on PC
1217f442ab6aSBarry Smith 
1218f442ab6aSBarry Smith    Input Parameters:
1219f442ab6aSBarry Smith .  pc - the preconditioner context
1220f442ab6aSBarry Smith 
1221f442ab6aSBarry Smith    Options Database Key:
1222f442ab6aSBarry Smith .  -pc_mg_distinct_smoothup
1223f442ab6aSBarry Smith 
1224f442ab6aSBarry Smith    Level: advanced
1225f442ab6aSBarry Smith 
122695452b02SPatrick Sanan    Notes:
122795452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
1228f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1229f442ab6aSBarry Smith 
1230f442ab6aSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1231f442ab6aSBarry Smith 
1232710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth()
1233f442ab6aSBarry Smith @*/
1234f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1235f442ab6aSBarry Smith {
1236f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1237f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1238f442ab6aSBarry Smith   PetscErrorCode ierr;
1239f442ab6aSBarry Smith   PetscInt       i,levels;
1240f442ab6aSBarry Smith   KSP            subksp;
1241f442ab6aSBarry Smith 
1242f442ab6aSBarry Smith   PetscFunctionBegin;
1243f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1244f442ab6aSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1245f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1246f442ab6aSBarry Smith 
1247f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1248710315b6SLawrence Mitchell     const char *prefix = NULL;
1249f442ab6aSBarry Smith     /* make sure smoother up and down are different */
1250f442ab6aSBarry Smith     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1251710315b6SLawrence Mitchell     ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr);
1252710315b6SLawrence Mitchell     ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr);
1253f442ab6aSBarry Smith     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1254f442ab6aSBarry Smith   }
1255f442ab6aSBarry Smith   PetscFunctionReturn(0);
1256f442ab6aSBarry Smith }
1257f442ab6aSBarry Smith 
12584b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
12594b9ad928SBarry Smith 
12603b09bd56SBarry Smith /*MC
1261ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
12623b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
12633b09bd56SBarry Smith 
12643b09bd56SBarry Smith    Options Database Keys:
12653b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1266391689abSStefano Zampini .  -pc_mg_cycle_type <v,w> - provide the cycle desired
12678c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
12683b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1269710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
12702134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
12718cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
12728cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1273e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
12748cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
12758cf6037eSBarry Smith                         to the binary output file called binaryoutput
12763b09bd56SBarry Smith 
127795452b02SPatrick Sanan    Notes:
12780b3c753eSRichard 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
12793b09bd56SBarry Smith 
12808cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
12818cf6037eSBarry Smith 
128223067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
128323067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
128423067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
128523067569SBarry Smith        residual is computed at the end of each cycle.
128623067569SBarry Smith 
12873b09bd56SBarry Smith    Level: intermediate
12883b09bd56SBarry Smith 
12898f87f92bSBarry Smith    Concepts: multigrid/multilevel
12903b09bd56SBarry Smith 
12918cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1292710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1293710315b6SLawrence Mitchell            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
129497177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
12950d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
12963b09bd56SBarry Smith M*/
12973b09bd56SBarry Smith 
12988cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
12994b9ad928SBarry Smith {
1300f3fbd535SBarry Smith   PC_MG          *mg;
1301f3fbd535SBarry Smith   PetscErrorCode ierr;
1302f3fbd535SBarry Smith 
13034b9ad928SBarry Smith   PetscFunctionBegin;
1304b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1305f3fbd535SBarry Smith   pc->data     = (void*)mg;
1306f3fbd535SBarry Smith   mg->nlevels  = -1;
130710eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
13082134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1309f3fbd535SBarry Smith 
131037a44384SMark Adams   pc->useAmat = PETSC_TRUE;
131137a44384SMark Adams 
13124b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
13134b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1314a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
13154b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
13164b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
13174b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1318fb15c04eSBarry Smith 
1319fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
132097d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
132197d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
13224b9ad928SBarry Smith   PetscFunctionReturn(0);
13234b9ad928SBarry Smith }
1324