xref: /petsc/src/ksp/pc/impls/mg/mg.c (revision d5a8f7e6b1c642fb14cbe2c9c1d7f6120d9aefd3)
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 
9f3b08a26SMatthew G. Knepley /*
10f3b08a26SMatthew G. Knepley    Contains the list of registered coarse space construction routines
11f3b08a26SMatthew G. Knepley */
12f3b08a26SMatthew G. Knepley PetscFunctionList PCMGCoarseList = NULL;
13f3b08a26SMatthew G. Knepley 
1431567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
154b9ad928SBarry Smith {
1631567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
1731567311SBarry Smith   PC_MG_Levels   *mgc,*mglevels = *mglevelsin;
186849ba73SBarry Smith   PetscErrorCode ierr;
1931567311SBarry Smith   PetscInt       cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
204b9ad928SBarry Smith 
214b9ad928SBarry Smith   PetscFunctionBegin;
2263e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2331567311SBarry Smith   ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr);  /* pre-smooth */
24c0decd05SBarry Smith   ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr);
2563e6d426SJed Brown   if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
2631567311SBarry Smith   if (mglevels->level) {  /* not the coarsest grid */
2763e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
2831567311SBarry Smith     ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr);
2963e6d426SJed Brown     if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);}
304b9ad928SBarry Smith 
314b9ad928SBarry Smith     /* if on finest level and have convergence criteria set */
3231567311SBarry Smith     if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
334b9ad928SBarry Smith       PetscReal rnorm;
3431567311SBarry Smith       ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr);
354b9ad928SBarry Smith       if (rnorm <= mg->ttol) {
3670441072SBarry Smith         if (rnorm < mg->abstol) {
374d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_ATOL;
3857622a8eSBarry Smith           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr);
394b9ad928SBarry Smith         } else {
404d0a8057SBarry Smith           *reason = PCRICHARDSON_CONVERGED_RTOL;
4157622a8eSBarry Smith           ierr    = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);CHKERRQ(ierr);
424b9ad928SBarry Smith         }
434b9ad928SBarry Smith         PetscFunctionReturn(0);
444b9ad928SBarry Smith       }
454b9ad928SBarry Smith     }
464b9ad928SBarry Smith 
4731567311SBarry Smith     mgc = *(mglevelsin - 1);
4863e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
4931567311SBarry Smith     ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr);
5063e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
51efb30889SBarry Smith     ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr);
524b9ad928SBarry Smith     while (cycles--) {
5331567311SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr);
544b9ad928SBarry Smith     }
5563e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5631567311SBarry Smith     ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr);
5763e6d426SJed Brown     if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5863e6d426SJed Brown     if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5931567311SBarry Smith     ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr);    /* post smooth */
60c0decd05SBarry Smith     ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr);
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;
146f3b08a26SMatthew G. Knepley   PetscInt       i,c,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++) {
165f3b08a26SMatthew G. Knepley       if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {ierr = VecDestroy(&mglevels[i]->coarseSpace[c]);CHKERRQ(ierr);}
166f3b08a26SMatthew G. Knepley       ierr = PetscFree(mglevels[i]->coarseSpace);CHKERRQ(ierr);
167f3b08a26SMatthew G. Knepley       mglevels[i]->coarseSpace = NULL;
1683aeaf226SBarry Smith       ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr);
1693aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1703aeaf226SBarry Smith         ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr);
1713aeaf226SBarry Smith       }
1723aeaf226SBarry Smith       ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr);
1733aeaf226SBarry Smith     }
174f3b08a26SMatthew G. Knepley     mg->Nc = 0;
1753aeaf226SBarry Smith   }
1763aeaf226SBarry Smith   PetscFunctionReturn(0);
1773aeaf226SBarry Smith }
1783aeaf226SBarry Smith 
17997d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
1804b9ad928SBarry Smith {
181dfbe8321SBarry Smith   PetscErrorCode ierr;
182f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
183ce94432eSBarry Smith   MPI_Comm       comm;
1843aeaf226SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
18510eca3edSLisandro Dalcin   PCMGType       mgtype     = mg->am;
18610167fecSBarry Smith   PetscInt       mgctype    = (PetscInt) PC_MG_CYCLE_V;
187f3fbd535SBarry Smith   PetscInt       i;
188f3fbd535SBarry Smith   PetscMPIInt    size;
189f3fbd535SBarry Smith   const char     *prefix;
190f3fbd535SBarry Smith   PC             ipc;
1913aeaf226SBarry Smith   PetscInt       n;
1924b9ad928SBarry Smith 
1934b9ad928SBarry Smith   PetscFunctionBegin;
1940700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
195548767e0SJed Brown   PetscValidLogicalCollectiveInt(pc,levels,2);
196548767e0SJed Brown   if (mg->nlevels == levels) PetscFunctionReturn(0);
1971a97d7b6SStefano Zampini   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
1983aeaf226SBarry Smith   if (mglevels) {
19910eca3edSLisandro Dalcin     mgctype = mglevels[0]->cycles;
2003aeaf226SBarry Smith     /* changing the number of levels so free up the previous stuff */
2013aeaf226SBarry Smith     ierr = PCReset_MG(pc);CHKERRQ(ierr);
2023aeaf226SBarry Smith     n    = mglevels[0]->levels;
2033aeaf226SBarry Smith     for (i=0; i<n; i++) {
2043aeaf226SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
2053aeaf226SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
2063aeaf226SBarry Smith       }
2073aeaf226SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
2083aeaf226SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
2093aeaf226SBarry Smith     }
2103aeaf226SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
2113aeaf226SBarry Smith   }
212f3fbd535SBarry Smith 
213f3fbd535SBarry Smith   mg->nlevels = levels;
214f3fbd535SBarry Smith 
215785e854fSJed Brown   ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr);
2163bb1ff40SBarry Smith   ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr);
217f3fbd535SBarry Smith 
218f3fbd535SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
219f3fbd535SBarry Smith 
220a9db87a2SMatthew G Knepley   mg->stageApply = 0;
221f3fbd535SBarry Smith   for (i=0; i<levels; i++) {
222b00a9115SJed Brown     ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr);
2232fa5cd67SKarl Rupp 
224f3fbd535SBarry Smith     mglevels[i]->level               = i;
225f3fbd535SBarry Smith     mglevels[i]->levels              = levels;
22610eca3edSLisandro Dalcin     mglevels[i]->cycles              = mgctype;
227336babb1SJed Brown     mg->default_smoothu              = 2;
228336babb1SJed Brown     mg->default_smoothd              = 2;
22963e6d426SJed Brown     mglevels[i]->eventsmoothsetup    = 0;
23063e6d426SJed Brown     mglevels[i]->eventsmoothsolve    = 0;
23163e6d426SJed Brown     mglevels[i]->eventresidual       = 0;
23263e6d426SJed Brown     mglevels[i]->eventinterprestrict = 0;
233f3fbd535SBarry Smith 
234f3fbd535SBarry Smith     if (comms) comm = comms[i];
235*d5a8f7e6SBarry Smith     if (comm != MPI_COMM_NULL) {
236f3fbd535SBarry Smith       ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr);
237422a814eSBarry Smith       ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr);
2380ee9a668SBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr);
2390ee9a668SBarry Smith       ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr);
2400ee9a668SBarry Smith       ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr);
2410ee9a668SBarry Smith       if (i || levels == 1) {
2420ee9a668SBarry Smith         char tprefix[128];
2430ee9a668SBarry Smith 
244336babb1SJed Brown         ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr);
2450059c7bdSJed Brown         ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr);
246336babb1SJed Brown         ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr);
247336babb1SJed Brown         ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr);
248336babb1SJed Brown         ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr);
2490ee9a668SBarry Smith         ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr);
250f3fbd535SBarry Smith 
251*d5a8f7e6SBarry Smith         ierr = PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i);CHKERRQ(ierr);
2520ee9a668SBarry Smith         ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr);
2530ee9a668SBarry Smith       } else {
254f3fbd535SBarry Smith         ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr);
255f3fbd535SBarry Smith 
2569dbfc187SHong Zhang         /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
257f3fbd535SBarry Smith         ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr);
258f3fbd535SBarry Smith         ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr);
259f3fbd535SBarry Smith         ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
260f3fbd535SBarry Smith         if (size > 1) {
261f3fbd535SBarry Smith           ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr);
262f3fbd535SBarry Smith         } else {
263f3fbd535SBarry Smith           ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr);
264f3fbd535SBarry Smith         }
265753b7fb9SBarry Smith         ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr);
266f3fbd535SBarry Smith       }
2673bb1ff40SBarry Smith       ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr);
268*d5a8f7e6SBarry Smith     }
269f3fbd535SBarry Smith     mglevels[i]->smoothu = mglevels[i]->smoothd;
27031567311SBarry Smith     mg->rtol             = 0.0;
27131567311SBarry Smith     mg->abstol           = 0.0;
27231567311SBarry Smith     mg->dtol             = 0.0;
27331567311SBarry Smith     mg->ttol             = 0.0;
27431567311SBarry Smith     mg->cyclesperpcapply = 1;
275f3fbd535SBarry Smith   }
276f3fbd535SBarry Smith   mg->levels = mglevels;
27710eca3edSLisandro Dalcin   ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
2784b9ad928SBarry Smith   PetscFunctionReturn(0);
2794b9ad928SBarry Smith }
2804b9ad928SBarry Smith 
28197d33e41SMatthew G. Knepley /*@C
28297d33e41SMatthew G. Knepley    PCMGSetLevels - Sets the number of levels to use with MG.
28397d33e41SMatthew G. Knepley    Must be called before any other MG routine.
28497d33e41SMatthew G. Knepley 
28597d33e41SMatthew G. Knepley    Logically Collective on PC
28697d33e41SMatthew G. Knepley 
28797d33e41SMatthew G. Knepley    Input Parameters:
28897d33e41SMatthew G. Knepley +  pc - the preconditioner context
28997d33e41SMatthew G. Knepley .  levels - the number of levels
29097d33e41SMatthew G. Knepley -  comms - optional communicators for each level; this is to allow solving the coarser problems
291*d5a8f7e6SBarry Smith            on smaller sets of processes. For processes that are not included in the computation
292*d5a8f7e6SBarry Smith            you must pass MPI_COMM_NULL.
29397d33e41SMatthew G. Knepley 
29497d33e41SMatthew G. Knepley    Level: intermediate
29597d33e41SMatthew G. Knepley 
29697d33e41SMatthew G. Knepley    Notes:
29797d33e41SMatthew G. Knepley      If the number of levels is one then the multigrid uses the -mg_levels prefix
29897d33e41SMatthew G. Knepley      for setting the level options rather than the -mg_coarse prefix.
29997d33e41SMatthew G. Knepley 
300*d5a8f7e6SBarry Smith      You can free the information in comms after this routine is called.
301*d5a8f7e6SBarry Smith 
302*d5a8f7e6SBarry Smith      The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level
303*d5a8f7e6SBarry Smith      are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on
304*d5a8f7e6SBarry Smith      the two levels, rank 0 in the original communicator will pass in an array of 2 communicators
305*d5a8f7e6SBarry Smith      of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators
306*d5a8f7e6SBarry Smith      the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate
307*d5a8f7e6SBarry Smith      in the coarse grid solve.
308*d5a8f7e6SBarry Smith 
309*d5a8f7e6SBarry Smith      Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one
310*d5a8f7e6SBarry Smith      must take special care in providing the restriction and interpolation operation. We recommend
311*d5a8f7e6SBarry Smith      providing these as two step operations; first perform a standard restriction or interpolation on
312*d5a8f7e6SBarry Smith      the full number of ranks for that level and then use an MPI call to copy the resulting vector
313*d5a8f7e6SBarry Smith      array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, not in both
314*d5a8f7e6SBarry Smith      cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and
315*d5a8f7e6SBarry Smith      recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries.
316*d5a8f7e6SBarry Smith 
317*d5a8f7e6SBarry Smith 
31897d33e41SMatthew G. Knepley .seealso: PCMGSetType(), PCMGGetLevels()
31997d33e41SMatthew G. Knepley @*/
32097d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
32197d33e41SMatthew G. Knepley {
32297d33e41SMatthew G. Knepley   PetscErrorCode ierr;
32397d33e41SMatthew G. Knepley 
32497d33e41SMatthew G. Knepley   PetscFunctionBegin;
32597d33e41SMatthew G. Knepley   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
32697d33e41SMatthew G. Knepley   if (comms) PetscValidPointer(comms,3);
32737b5128cSMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr);
32897d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
32997d33e41SMatthew G. Knepley }
33097d33e41SMatthew G. Knepley 
331c07bf074SBarry Smith 
332c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc)
333c07bf074SBarry Smith {
334c07bf074SBarry Smith   PetscErrorCode ierr;
335a06653b4SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
336a06653b4SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
337a06653b4SBarry Smith   PetscInt       i,n;
338c07bf074SBarry Smith 
339c07bf074SBarry Smith   PetscFunctionBegin;
340a06653b4SBarry Smith   ierr = PCReset_MG(pc);CHKERRQ(ierr);
341a06653b4SBarry Smith   if (mglevels) {
342a06653b4SBarry Smith     n = mglevels[0]->levels;
343a06653b4SBarry Smith     for (i=0; i<n; i++) {
344a06653b4SBarry Smith       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
3456bf464f9SBarry Smith         ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr);
346a06653b4SBarry Smith       }
3476bf464f9SBarry Smith       ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr);
348a06653b4SBarry Smith       ierr = PetscFree(mglevels[i]);CHKERRQ(ierr);
349a06653b4SBarry Smith     }
350a06653b4SBarry Smith     ierr = PetscFree(mg->levels);CHKERRQ(ierr);
351a06653b4SBarry Smith   }
352c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
353fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);CHKERRQ(ierr);
354fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);CHKERRQ(ierr);
355f3fbd535SBarry Smith   PetscFunctionReturn(0);
356f3fbd535SBarry Smith }
357f3fbd535SBarry Smith 
358f3fbd535SBarry Smith 
359f3fbd535SBarry Smith 
36009573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
36109573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
36209573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
363f3fbd535SBarry Smith 
364f3fbd535SBarry Smith /*
365f3fbd535SBarry Smith    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
366f3fbd535SBarry Smith              or full cycle of multigrid.
367f3fbd535SBarry Smith 
368f3fbd535SBarry Smith   Note:
369f3fbd535SBarry Smith   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
370f3fbd535SBarry Smith */
371f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
372f3fbd535SBarry Smith {
373f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
374f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
375f3fbd535SBarry Smith   PetscErrorCode ierr;
376391689abSStefano Zampini   PC             tpc;
377f3fbd535SBarry Smith   PetscInt       levels = mglevels[0]->levels,i;
378391689abSStefano Zampini   PetscBool      changeu,changed;
379f3fbd535SBarry Smith 
380f3fbd535SBarry Smith   PetscFunctionBegin;
381a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);}
382e1d8e5deSBarry Smith   /* When the DM is supplying the matrix then it will not exist until here */
383298cc208SBarry Smith   for (i=0; i<levels; i++) {
384e1d8e5deSBarry Smith     if (!mglevels[i]->A) {
38523ee1639SBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr);
386298cc208SBarry Smith       ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr);
387e1d8e5deSBarry Smith     }
388e1d8e5deSBarry Smith   }
389e1d8e5deSBarry Smith 
390391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr);
391391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr);
392391689abSStefano Zampini   ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr);
393391689abSStefano Zampini   ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr);
394391689abSStefano Zampini   if (!changeu && !changed) {
395391689abSStefano Zampini     ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr);
396f3fbd535SBarry Smith     mglevels[levels-1]->b = b;
397391689abSStefano Zampini   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
398391689abSStefano Zampini     if (!mglevels[levels-1]->b) {
399391689abSStefano Zampini       Vec *vec;
400391689abSStefano Zampini 
401391689abSStefano Zampini       ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
402391689abSStefano Zampini       mglevels[levels-1]->b = *vec;
4037ae21d43SStefano Zampini       ierr = PetscFree(vec);CHKERRQ(ierr);
404391689abSStefano Zampini     }
405391689abSStefano Zampini     ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr);
406391689abSStefano Zampini   }
407f3fbd535SBarry Smith   mglevels[levels-1]->x = x;
408391689abSStefano Zampini 
40931567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
410f3fbd535SBarry Smith     ierr = VecSet(x,0.0);CHKERRQ(ierr);
41131567311SBarry Smith     for (i=0; i<mg->cyclesperpcapply; i++) {
4120298fd71SBarry Smith       ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr);
413f3fbd535SBarry Smith     }
4142fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_ADDITIVE) {
41531567311SBarry Smith     ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr);
4162fa5cd67SKarl Rupp   } else if (mg->am == PC_MG_KASKADE) {
41731567311SBarry Smith     ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr);
4182fa5cd67SKarl Rupp   } else {
419f3fbd535SBarry Smith     ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr);
420f3fbd535SBarry Smith   }
421a9db87a2SMatthew G Knepley   if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);}
422391689abSStefano Zampini   if (!changeu && !changed) mglevels[levels-1]->b = NULL;
423f3fbd535SBarry Smith   PetscFunctionReturn(0);
424f3fbd535SBarry Smith }
425f3fbd535SBarry Smith 
426f3fbd535SBarry Smith 
4274416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
428f3fbd535SBarry Smith {
429f3fbd535SBarry Smith   PetscErrorCode   ierr;
430710315b6SLawrence Mitchell   PetscInt         levels,cycles;
431f3b08a26SMatthew G. Knepley   PetscBool        flg, flg2;
432f3fbd535SBarry Smith   PC_MG            *mg = (PC_MG*)pc->data;
4333d3eaba7SBarry Smith   PC_MG_Levels     **mglevels;
434f3fbd535SBarry Smith   PCMGType         mgtype;
435f3fbd535SBarry Smith   PCMGCycleType    mgctype;
4362134b1e4SBarry Smith   PCMGGalerkinType gtype;
437f3fbd535SBarry Smith 
438f3fbd535SBarry Smith   PetscFunctionBegin;
4391d743356SStefano Zampini   levels = PetscMax(mg->nlevels,1);
440e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr);
441f3fbd535SBarry Smith   ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr);
4421a97d7b6SStefano Zampini   if (!flg && !mg->levels && pc->dm) {
443eb3f98d2SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
444eb3f98d2SBarry Smith     levels++;
4453aeaf226SBarry Smith     mg->usedmfornumberoflevels = PETSC_TRUE;
446eb3f98d2SBarry Smith   }
4470298fd71SBarry Smith   ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
4483aeaf226SBarry Smith   mglevels = mg->levels;
4493aeaf226SBarry Smith 
450f3fbd535SBarry Smith   mgctype = (PCMGCycleType) mglevels[0]->cycles;
451f3fbd535SBarry Smith   ierr    = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr);
452f3fbd535SBarry Smith   if (flg) {
453f3fbd535SBarry Smith     ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr);
4542fa5cd67SKarl Rupp   }
4552134b1e4SBarry Smith   gtype = mg->galerkin;
4562134b1e4SBarry Smith   ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)&gtype,&flg);CHKERRQ(ierr);
4572134b1e4SBarry Smith   if (flg) {
4582134b1e4SBarry Smith     ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr);
459f3fbd535SBarry Smith   }
460f3b08a26SMatthew G. Knepley   flg2 = PETSC_FALSE;
461f3b08a26SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr);
462f3b08a26SMatthew G. Knepley   if (flg) {ierr = PCMGSetAdaptInterpolation(pc, flg2);CHKERRQ(ierr);}
463f3b08a26SMatthew G. Knepley   ierr = PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg);CHKERRQ(ierr);
464f3b08a26SMatthew G. Knepley   ierr = PetscOptionsEnum("-pc_mg_adapt_interp_coarse_space","Type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector","PCMGSetAdaptCoarseSpaceType",PCMGCoarseSpaceTypes,(PetscEnum)mg->coarseSpaceType,(PetscEnum*)&mg->coarseSpaceType,&flg);CHKERRQ(ierr);
465f3b08a26SMatthew G. Knepley   ierr = PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);CHKERRQ(ierr);
466f56b1265SBarry Smith   flg = PETSC_FALSE;
4678e5aa403SBarry Smith   ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr);
468f442ab6aSBarry Smith   if (flg) {
469f442ab6aSBarry Smith     ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr);
470f442ab6aSBarry Smith   }
47131567311SBarry Smith   mgtype = mg->am;
472f3fbd535SBarry Smith   ierr   = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr);
473f3fbd535SBarry Smith   if (flg) {
474f3fbd535SBarry Smith     ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr);
475f3fbd535SBarry Smith   }
47631567311SBarry Smith   if (mg->am == PC_MG_MULTIPLICATIVE) {
4775363c1e0SLisandro Dalcin     ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr);
478f3fbd535SBarry Smith     if (flg) {
479f3fbd535SBarry Smith       ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr);
480f3fbd535SBarry Smith     }
481f3fbd535SBarry Smith   }
482f3fbd535SBarry Smith   flg  = PETSC_FALSE;
4830298fd71SBarry Smith   ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr);
484f3fbd535SBarry Smith   if (flg) {
485f3fbd535SBarry Smith     PetscInt i;
486f3fbd535SBarry Smith     char     eventname[128];
4871a97d7b6SStefano Zampini 
488f3fbd535SBarry Smith     levels = mglevels[0]->levels;
489f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
490f3fbd535SBarry Smith       sprintf(eventname,"MGSetup Level %d",(int)i);
49163e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr);
492f3fbd535SBarry Smith       sprintf(eventname,"MGSmooth Level %d",(int)i);
49363e6d426SJed Brown       ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr);
494f3fbd535SBarry Smith       if (i) {
495f3fbd535SBarry Smith         sprintf(eventname,"MGResid Level %d",(int)i);
49663e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr);
497f3fbd535SBarry Smith         sprintf(eventname,"MGInterp Level %d",(int)i);
49863e6d426SJed Brown         ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr);
499f3fbd535SBarry Smith       }
500f3fbd535SBarry Smith     }
501eec35531SMatthew G Knepley 
502ec60ca73SSatish Balay #if defined(PETSC_USE_LOG)
503eec35531SMatthew G Knepley     {
504eec35531SMatthew G Knepley       const char    *sname = "MG Apply";
505eec35531SMatthew G Knepley       PetscStageLog stageLog;
506eec35531SMatthew G Knepley       PetscInt      st;
507eec35531SMatthew G Knepley 
508eec35531SMatthew G Knepley       ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr);
509eec35531SMatthew G Knepley       for (st = 0; st < stageLog->numStages; ++st) {
510eec35531SMatthew G Knepley         PetscBool same;
511eec35531SMatthew G Knepley 
512eec35531SMatthew G Knepley         ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr);
5132fa5cd67SKarl Rupp         if (same) mg->stageApply = st;
514eec35531SMatthew G Knepley       }
515eec35531SMatthew G Knepley       if (!mg->stageApply) {
516eec35531SMatthew G Knepley         ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr);
517eec35531SMatthew G Knepley       }
518eec35531SMatthew G Knepley     }
519ec60ca73SSatish Balay #endif
520f3fbd535SBarry Smith   }
521f3fbd535SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
522f3b08a26SMatthew G. Knepley   /* Check option consistency */
523f3b08a26SMatthew G. Knepley   ierr = PCMGGetGalerkin(pc, &gtype);CHKERRQ(ierr);
524f3b08a26SMatthew G. Knepley   ierr = PCMGGetAdaptInterpolation(pc, &flg);CHKERRQ(ierr);
525f3b08a26SMatthew G. Knepley   if (flg && (gtype >= PC_MG_GALERKIN_NONE)) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
526f3fbd535SBarry Smith   PetscFunctionReturn(0);
527f3fbd535SBarry Smith }
528f3fbd535SBarry Smith 
5290a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL};
5300a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL};
5310a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL};
532f3b08a26SMatthew G. Knepley const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL};
533f3fbd535SBarry Smith 
5349804daf3SBarry Smith #include <petscdraw.h>
535f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
536f3fbd535SBarry Smith {
537f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
538f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
539f3fbd535SBarry Smith   PetscErrorCode ierr;
540e3deeeafSJed Brown   PetscInt       levels = mglevels ? mglevels[0]->levels : 0,i;
541219b1045SBarry Smith   PetscBool      iascii,isbinary,isdraw;
542f3fbd535SBarry Smith 
543f3fbd535SBarry Smith   PetscFunctionBegin;
544251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
5455b0b0462SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
546219b1045SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
547f3fbd535SBarry Smith   if (iascii) {
548e3deeeafSJed Brown     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
549efd4aadfSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr);
55031567311SBarry Smith     if (mg->am == PC_MG_MULTIPLICATIVE) {
55131567311SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr);
552f3fbd535SBarry Smith     }
5532134b1e4SBarry Smith     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
554f3fbd535SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
5552134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
5562134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr);
5572134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
5582134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr);
5592134b1e4SBarry Smith     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
5602134b1e4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr);
5614f66f45eSBarry Smith     } else {
5624f66f45eSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr);
563f3fbd535SBarry Smith     }
5645adeb434SBarry Smith     if (mg->view){
5655adeb434SBarry Smith       ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr);
5665adeb434SBarry Smith     }
567f3fbd535SBarry Smith     for (i=0; i<levels; i++) {
568f3fbd535SBarry Smith       if (!i) {
56989cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr);
570f3fbd535SBarry Smith       } else {
57189cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
572f3fbd535SBarry Smith       }
573f3fbd535SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
574f3fbd535SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
575f3fbd535SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
576f3fbd535SBarry Smith       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
577f3fbd535SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
578f3fbd535SBarry Smith       } else if (i) {
57989cce641SBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
580f3fbd535SBarry Smith         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
581f3fbd535SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
582f3fbd535SBarry Smith         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
583f3fbd535SBarry Smith       }
584f3fbd535SBarry Smith     }
5855b0b0462SBarry Smith   } else if (isbinary) {
5865b0b0462SBarry Smith     for (i=levels-1; i>=0; i--) {
5875b0b0462SBarry Smith       ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
5885b0b0462SBarry Smith       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
5895b0b0462SBarry Smith         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
5905b0b0462SBarry Smith       }
5915b0b0462SBarry Smith     }
592219b1045SBarry Smith   } else if (isdraw) {
593219b1045SBarry Smith     PetscDraw draw;
594b4375e8dSPeter Brune     PetscReal x,w,y,bottom,th;
595219b1045SBarry Smith     ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
596219b1045SBarry Smith     ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
5970298fd71SBarry Smith     ierr   = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr);
598219b1045SBarry Smith     bottom = y - th;
599219b1045SBarry Smith     for (i=levels-1; i>=0; i--) {
600b4375e8dSPeter Brune       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
601219b1045SBarry Smith         ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
602219b1045SBarry Smith         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
603219b1045SBarry Smith         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
604b4375e8dSPeter Brune       } else {
605b4375e8dSPeter Brune         w    = 0.5*PetscMin(1.0-x,x);
606b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
607b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr);
608b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
609b4375e8dSPeter Brune         ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
610b4375e8dSPeter Brune         ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr);
611b4375e8dSPeter Brune         ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
612b4375e8dSPeter Brune       }
6130298fd71SBarry Smith       ierr    = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr);
6141cd381d6SBarry Smith       bottom -= th;
615219b1045SBarry Smith     }
6165b0b0462SBarry Smith   }
617f3fbd535SBarry Smith   PetscFunctionReturn(0);
618f3fbd535SBarry Smith }
619f3fbd535SBarry Smith 
620af0996ceSBarry Smith #include <petsc/private/kspimpl.h>
621cab2e9ccSBarry Smith 
622f3fbd535SBarry Smith /*
623f3fbd535SBarry Smith     Calls setup for the KSP on each level
624f3fbd535SBarry Smith */
625f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc)
626f3fbd535SBarry Smith {
627f3fbd535SBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
628f3fbd535SBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
629f3fbd535SBarry Smith   PetscErrorCode ierr;
6301a97d7b6SStefano Zampini   PetscInt       i,n;
63198e04984SBarry Smith   PC             cpc;
6323db492dfSBarry Smith   PetscBool      dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
633f3fbd535SBarry Smith   Mat            dA,dB;
634f3fbd535SBarry Smith   Vec            tvec;
635218a07d4SBarry Smith   DM             *dms;
6360a545947SLisandro Dalcin   PetscViewer    viewer = NULL;
6372134b1e4SBarry Smith   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE;
638f3fbd535SBarry Smith 
639f3fbd535SBarry Smith   PetscFunctionBegin;
6401a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
6411a97d7b6SStefano Zampini   n = mglevels[0]->levels;
64201bc414fSDmitry Karpeev   /* FIX: Move this to PCSetFromOptions_MG? */
6433aeaf226SBarry Smith   if (mg->usedmfornumberoflevels) {
6443aeaf226SBarry Smith     PetscInt levels;
6453aeaf226SBarry Smith     ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr);
6463aeaf226SBarry Smith     levels++;
6473aeaf226SBarry Smith     if (levels > n) { /* the problem is now being solved on a finer grid */
6480298fd71SBarry Smith       ierr     = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr);
6493aeaf226SBarry Smith       n        = levels;
6503aeaf226SBarry Smith       ierr     = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
6513aeaf226SBarry Smith       mglevels = mg->levels;
6523aeaf226SBarry Smith     }
6533aeaf226SBarry Smith   }
65498e04984SBarry Smith   ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr);
6553aeaf226SBarry Smith 
656f3fbd535SBarry Smith 
657f3fbd535SBarry Smith   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
658f3fbd535SBarry Smith   /* so use those from global PC */
659f3fbd535SBarry Smith   /* Is this what we always want? What if user wants to keep old one? */
6600298fd71SBarry Smith   ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr);
66198e04984SBarry Smith   if (opsset) {
66298e04984SBarry Smith     Mat mmat;
66323ee1639SBarry Smith     ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr);
66498e04984SBarry Smith     if (mmat == pc->pmat) opsset = PETSC_FALSE;
66598e04984SBarry Smith   }
6664a5f07a5SMark F. Adams 
6674a5f07a5SMark F. Adams   if (!opsset) {
66871b23a65SMark F. Adams     ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr);
6694a5f07a5SMark F. Adams     if (use_amat) {
670f3fbd535SBarry 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);
67123ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr);
67269858f1bSStefano Zampini     } else {
6734a5f07a5SMark 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);
67423ee1639SBarry Smith       ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr);
6754a5f07a5SMark F. Adams     }
6764a5f07a5SMark F. Adams   }
677f3fbd535SBarry Smith 
6789c7ffb39SBarry Smith   for (i=n-1; i>0; i--) {
6799c7ffb39SBarry Smith     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
6809c7ffb39SBarry Smith       missinginterpolate = PETSC_TRUE;
6819c7ffb39SBarry Smith       continue;
6829c7ffb39SBarry Smith     }
6839c7ffb39SBarry Smith   }
6842134b1e4SBarry Smith 
6852134b1e4SBarry Smith   ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr);
6862134b1e4SBarry Smith   if (dA == dB) dAeqdB = PETSC_TRUE;
6872134b1e4SBarry Smith   if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
6882134b1e4SBarry Smith     needRestricts = PETSC_TRUE;  /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
6892134b1e4SBarry Smith   }
6902134b1e4SBarry Smith 
6912134b1e4SBarry Smith 
6929c7ffb39SBarry Smith   /*
6939c7ffb39SBarry Smith    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
6949c7ffb39SBarry Smith    Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
6959c7ffb39SBarry Smith   */
6962134b1e4SBarry Smith   if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
6972d2b81a6SBarry Smith         /* construct the interpolation from the DMs */
6982e499ae9SBarry Smith     Mat p;
69973250ac0SBarry Smith     Vec rscale;
700785e854fSJed Brown     ierr     = PetscMalloc1(n,&dms);CHKERRQ(ierr);
701218a07d4SBarry Smith     dms[n-1] = pc->dm;
702ef1267afSMatthew G. Knepley     /* Separately create them so we do not get DMKSP interference between levels */
703ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);}
70441fce8fdSHong Zhang         /*
70541fce8fdSHong Zhang            Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
70641fce8fdSHong Zhang            Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
70741fce8fdSHong Zhang            But it is safe to use -dm_mat_type.
70841fce8fdSHong Zhang 
70941fce8fdSHong Zhang            The mat type should not be hardcoded like this, we need to find a better way.
71041fce8fdSHong Zhang     ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr);
71141fce8fdSHong Zhang     */
712218a07d4SBarry Smith     for (i=n-2; i>-1; i--) {
713942e3340SBarry Smith       DMKSP     kdm;
714eab52d2bSLawrence Mitchell       PetscBool dmhasrestrict, dmhasinject;
7153c0c59f3SBarry Smith       ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr);
7162134b1e4SBarry Smith       if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);}
717c27ee7a3SPatrick Farrell       if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
718c27ee7a3SPatrick Farrell         ierr = KSPSetDM(mglevels[i]->smoothu,dms[i]);CHKERRQ(ierr);
719c27ee7a3SPatrick Farrell         if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);CHKERRQ(ierr);}
720c27ee7a3SPatrick Farrell       }
721942e3340SBarry Smith       ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr);
722d1a3e677SJed Brown       /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
723d1a3e677SJed Brown        * a bitwise OR of computing the matrix, RHS, and initial iterate. */
7240298fd71SBarry Smith       kdm->ops->computerhs = NULL;
7250298fd71SBarry Smith       kdm->rhsctx          = NULL;
72624c3aa18SBarry Smith       if (!mglevels[i+1]->interpolate) {
727e727c939SJed Brown         ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr);
7282d2b81a6SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr);
72900ac8be1SBarry Smith         if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);}
73073250ac0SBarry Smith         ierr = VecDestroy(&rscale);CHKERRQ(ierr);
7316bf464f9SBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
732218a07d4SBarry Smith       }
7333ad4599aSBarry Smith       ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr);
7343ad4599aSBarry Smith       if (dmhasrestrict && !mglevels[i+1]->restrct){
7353ad4599aSBarry Smith         ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr);
7363ad4599aSBarry Smith         ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr);
7373ad4599aSBarry Smith         ierr = MatDestroy(&p);CHKERRQ(ierr);
7383ad4599aSBarry Smith       }
739eab52d2bSLawrence Mitchell       ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr);
740eab52d2bSLawrence Mitchell       if (dmhasinject && !mglevels[i+1]->inject){
741eab52d2bSLawrence Mitchell         ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr);
742eab52d2bSLawrence Mitchell         ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr);
743eab52d2bSLawrence Mitchell         ierr = MatDestroy(&p);CHKERRQ(ierr);
744eab52d2bSLawrence Mitchell       }
74524c3aa18SBarry Smith     }
7462d2b81a6SBarry Smith 
747ef1267afSMatthew G. Knepley     for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);}
7482d2b81a6SBarry Smith     ierr = PetscFree(dms);CHKERRQ(ierr);
749ce4cda84SJed Brown   }
750cab2e9ccSBarry Smith 
751ce4cda84SJed Brown   if (pc->dm && !pc->setupcalled) {
7522134b1e4SBarry Smith     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
753cab2e9ccSBarry Smith     ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr);
754cab2e9ccSBarry Smith     ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr);
755c27ee7a3SPatrick Farrell     if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) {
756c27ee7a3SPatrick Farrell       ierr = KSPSetDM(mglevels[n-1]->smoothu,pc->dm);CHKERRQ(ierr);
757c27ee7a3SPatrick Farrell       ierr = KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);CHKERRQ(ierr);
758c27ee7a3SPatrick Farrell     }
759218a07d4SBarry Smith   }
760218a07d4SBarry Smith 
7612134b1e4SBarry Smith   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
7622134b1e4SBarry Smith     Mat       A,B;
7632134b1e4SBarry Smith     PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
7642134b1e4SBarry Smith     MatReuse  reuse = MAT_INITIAL_MATRIX;
7652134b1e4SBarry Smith 
7662134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
7672134b1e4SBarry Smith     if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
7682134b1e4SBarry Smith     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
769f3fbd535SBarry Smith     for (i=n-2; i>-1; i--) {
770b9367914SBarry 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");
771b9367914SBarry Smith       if (!mglevels[i+1]->interpolate) {
772b9367914SBarry Smith         ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr);
773b9367914SBarry Smith       }
774b9367914SBarry Smith       if (!mglevels[i+1]->restrct) {
775b9367914SBarry Smith         ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr);
776b9367914SBarry Smith       }
7772134b1e4SBarry Smith       if (reuse == MAT_REUSE_MATRIX) {
7782134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr);
7792134b1e4SBarry Smith       }
7802134b1e4SBarry Smith       if (doA) {
7812df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr);
7822134b1e4SBarry Smith       }
7832134b1e4SBarry Smith       if (doB) {
7842df6c5c3SPatrick Farrell         ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr);
785b9367914SBarry Smith       }
7862134b1e4SBarry Smith       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
7872134b1e4SBarry Smith       if (!doA && dAeqdB) {
7882134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);}
7892134b1e4SBarry Smith         A = B;
7902134b1e4SBarry Smith       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
7912134b1e4SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr);
7922134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
793b9367914SBarry Smith       }
7942134b1e4SBarry Smith       if (!doB && dAeqdB) {
7952134b1e4SBarry Smith         if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);}
7962134b1e4SBarry Smith         B = A;
7972134b1e4SBarry Smith       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
79823ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr);
7992134b1e4SBarry Smith         ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
8007d537d92SJed Brown       }
8012134b1e4SBarry Smith       if (reuse == MAT_INITIAL_MATRIX) {
8022134b1e4SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr);
8032134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr);
8042134b1e4SBarry Smith         ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr);
8052134b1e4SBarry Smith       }
8062134b1e4SBarry Smith       dA = A;
807f3fbd535SBarry Smith       dB = B;
808f3fbd535SBarry Smith     }
809f3fbd535SBarry Smith   }
810f3b08a26SMatthew G. Knepley 
811f3b08a26SMatthew G. Knepley 
812f3b08a26SMatthew G. Knepley   /* Adapt interpolation matrices */
813f3b08a26SMatthew G. Knepley   if (mg->adaptInterpolation) {
814f3b08a26SMatthew G. Knepley     mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */
815f3b08a26SMatthew G. Knepley     for (i = 0; i < n; ++i) {
816f3b08a26SMatthew G. Knepley       ierr = PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);CHKERRQ(ierr);
817f3b08a26SMatthew G. Knepley       if (i) {ierr = PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mg->Nc, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace);CHKERRQ(ierr);}
818f3b08a26SMatthew G. Knepley     }
819f3b08a26SMatthew G. Knepley     for (i = n-2; i > -1; --i) {
820f3b08a26SMatthew G. Knepley       ierr = PCMGRecomputeLevelOperators_Internal(pc, i);CHKERRQ(ierr);
821f3b08a26SMatthew G. Knepley     }
822f3b08a26SMatthew G. Knepley   }
823f3b08a26SMatthew G. Knepley 
8242134b1e4SBarry Smith   if (needRestricts && pc->dm) {
825caa4e7f2SJed Brown     for (i=n-2; i>=0; i--) {
826ccceb50cSJed Brown       DM  dmfine,dmcoarse;
827caa4e7f2SJed Brown       Mat Restrict,Inject;
828caa4e7f2SJed Brown       Vec rscale;
829ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr);
830ccceb50cSJed Brown       ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr);
831caa4e7f2SJed Brown       ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr);
832caa4e7f2SJed Brown       ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr);
833eab52d2bSLawrence Mitchell       ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr);
834caa4e7f2SJed Brown       ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr);
835caa4e7f2SJed Brown     }
836caa4e7f2SJed Brown   }
837f3fbd535SBarry Smith 
838f3fbd535SBarry Smith   if (!pc->setupcalled) {
839f3fbd535SBarry Smith     for (i=0; i<n; i++) {
840f3fbd535SBarry Smith       ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr);
841f3fbd535SBarry Smith     }
842f3fbd535SBarry Smith     for (i=1; i<n; i++) {
843f3fbd535SBarry Smith       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
844f3fbd535SBarry Smith         ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr);
845f3fbd535SBarry Smith       }
846f3fbd535SBarry Smith     }
8473ad4599aSBarry Smith     /* insure that if either interpolation or restriction is set the other other one is set */
848f3fbd535SBarry Smith     for (i=1; i<n; i++) {
8493ad4599aSBarry Smith       ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr);
8503ad4599aSBarry Smith       ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr);
851f3fbd535SBarry Smith     }
852f3fbd535SBarry Smith     for (i=0; i<n-1; i++) {
853f3fbd535SBarry Smith       if (!mglevels[i]->b) {
854f3fbd535SBarry Smith         Vec *vec;
8552a7a6963SBarry Smith         ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
856f3fbd535SBarry Smith         ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr);
8576bf464f9SBarry Smith         ierr = VecDestroy(vec);CHKERRQ(ierr);
858f3fbd535SBarry Smith         ierr = PetscFree(vec);CHKERRQ(ierr);
859f3fbd535SBarry Smith       }
860f3fbd535SBarry Smith       if (!mglevels[i]->r && i) {
861f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
862f3fbd535SBarry Smith         ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr);
8636bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
864f3fbd535SBarry Smith       }
865f3fbd535SBarry Smith       if (!mglevels[i]->x) {
866f3fbd535SBarry Smith         ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr);
867f3fbd535SBarry Smith         ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr);
8686bf464f9SBarry Smith         ierr = VecDestroy(&tvec);CHKERRQ(ierr);
869f3fbd535SBarry Smith       }
870f3fbd535SBarry Smith     }
871f3fbd535SBarry Smith     if (n != 1 && !mglevels[n-1]->r) {
872f3fbd535SBarry Smith       /* PCMGSetR() on the finest level if user did not supply it */
873f3fbd535SBarry Smith       Vec *vec;
8742a7a6963SBarry Smith       ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr);
875f3fbd535SBarry Smith       ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr);
8766bf464f9SBarry Smith       ierr = VecDestroy(vec);CHKERRQ(ierr);
877f3fbd535SBarry Smith       ierr = PetscFree(vec);CHKERRQ(ierr);
878f3fbd535SBarry Smith     }
879f3fbd535SBarry Smith   }
880f3fbd535SBarry Smith 
88198e04984SBarry Smith   if (pc->dm) {
88298e04984SBarry Smith     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
88398e04984SBarry Smith     for (i=0; i<n-1; i++) {
88498e04984SBarry Smith       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
88598e04984SBarry Smith     }
88698e04984SBarry Smith   }
887f3fbd535SBarry Smith 
888f3fbd535SBarry Smith   for (i=1; i<n; i++) {
8892a21e185SBarry Smith     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
890f3fbd535SBarry Smith       /* if doing only down then initial guess is zero */
891f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr);
892f3fbd535SBarry Smith     }
89363e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
894f3fbd535SBarry Smith     ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr);
895c0decd05SBarry Smith     if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
896899639b0SHong Zhang       pc->failedreason = PC_SUBPC_ERROR;
897899639b0SHong Zhang     }
89863e6d426SJed Brown     if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
899d42688cbSBarry Smith     if (!mglevels[i]->residual) {
900d42688cbSBarry Smith       Mat mat;
90113842ffbSBarry Smith       ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr);
90254b2cd4bSJed Brown       ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr);
903d42688cbSBarry Smith     }
904f3fbd535SBarry Smith   }
905f3fbd535SBarry Smith   for (i=1; i<n; i++) {
906f3fbd535SBarry Smith     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
907f3fbd535SBarry Smith       Mat downmat,downpmat;
908f3fbd535SBarry Smith 
909f3fbd535SBarry Smith       /* check if operators have been set for up, if not use down operators to set them */
9100298fd71SBarry Smith       ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr);
911f3fbd535SBarry Smith       if (!opsset) {
91223ee1639SBarry Smith         ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr);
91323ee1639SBarry Smith         ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr);
914f3fbd535SBarry Smith       }
915f3fbd535SBarry Smith 
916f3fbd535SBarry Smith       ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr);
91763e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
918f3fbd535SBarry Smith       ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr);
919c0decd05SBarry Smith       if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
920899639b0SHong Zhang         pc->failedreason = PC_SUBPC_ERROR;
921899639b0SHong Zhang       }
92263e6d426SJed Brown       if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
923f3fbd535SBarry Smith     }
924f3fbd535SBarry Smith   }
925f3fbd535SBarry Smith 
92663e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
927f3fbd535SBarry Smith   ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr);
928c0decd05SBarry Smith   if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
929899639b0SHong Zhang     pc->failedreason = PC_SUBPC_ERROR;
930899639b0SHong Zhang   }
93163e6d426SJed Brown   if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
932f3fbd535SBarry Smith 
933f3fbd535SBarry Smith   /*
934f3fbd535SBarry Smith      Dump the interpolation/restriction matrices plus the
935e3c5b3baSBarry Smith    Jacobian/stiffness on each level. This allows MATLAB users to
936f3fbd535SBarry Smith    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
937f3fbd535SBarry Smith 
938f3fbd535SBarry Smith    Only support one or the other at the same time.
939f3fbd535SBarry Smith   */
940f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER)
941c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr);
942ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
943f3fbd535SBarry Smith   dump = PETSC_FALSE;
944f3fbd535SBarry Smith #endif
945c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr);
946ce94432eSBarry Smith   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
947f3fbd535SBarry Smith 
948f3fbd535SBarry Smith   if (viewer) {
949f3fbd535SBarry Smith     for (i=1; i<n; i++) {
950f3fbd535SBarry Smith       ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr);
951f3fbd535SBarry Smith     }
952f3fbd535SBarry Smith     for (i=0; i<n; i++) {
953f3fbd535SBarry Smith       ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr);
954f3fbd535SBarry Smith       ierr = MatView(pc->mat,viewer);CHKERRQ(ierr);
955f3fbd535SBarry Smith     }
956f3fbd535SBarry Smith   }
957f3fbd535SBarry Smith   PetscFunctionReturn(0);
958f3fbd535SBarry Smith }
959f3fbd535SBarry Smith 
960f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/
961f3fbd535SBarry Smith 
96297d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
96397d33e41SMatthew G. Knepley {
96497d33e41SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
96597d33e41SMatthew G. Knepley 
96697d33e41SMatthew G. Knepley   PetscFunctionBegin;
96797d33e41SMatthew G. Knepley   *levels = mg->nlevels;
96897d33e41SMatthew G. Knepley   PetscFunctionReturn(0);
96997d33e41SMatthew G. Knepley }
97097d33e41SMatthew G. Knepley 
9714b9ad928SBarry Smith /*@
97297177400SBarry Smith    PCMGGetLevels - Gets the number of levels to use with MG.
9734b9ad928SBarry Smith 
9744b9ad928SBarry Smith    Not Collective
9754b9ad928SBarry Smith 
9764b9ad928SBarry Smith    Input Parameter:
9774b9ad928SBarry Smith .  pc - the preconditioner context
9784b9ad928SBarry Smith 
9794b9ad928SBarry Smith    Output parameter:
9804b9ad928SBarry Smith .  levels - the number of levels
9814b9ad928SBarry Smith 
9824b9ad928SBarry Smith    Level: advanced
9834b9ad928SBarry Smith 
98497177400SBarry Smith .seealso: PCMGSetLevels()
9854b9ad928SBarry Smith @*/
9867087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
9874b9ad928SBarry Smith {
988e952c529SMatthew G. Knepley   PetscErrorCode ierr;
9894b9ad928SBarry Smith 
9904b9ad928SBarry Smith   PetscFunctionBegin;
9910700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
9924482741eSBarry Smith   PetscValidIntPointer(levels,2);
99397d33e41SMatthew G. Knepley   *levels = 0;
99437b5128cSMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr);
9954b9ad928SBarry Smith   PetscFunctionReturn(0);
9964b9ad928SBarry Smith }
9974b9ad928SBarry Smith 
9984b9ad928SBarry Smith /*@
99997177400SBarry Smith    PCMGSetType - Determines the form of multigrid to use:
10004b9ad928SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
10014b9ad928SBarry Smith 
1002ad4df100SBarry Smith    Logically Collective on PC
10034b9ad928SBarry Smith 
10044b9ad928SBarry Smith    Input Parameters:
10054b9ad928SBarry Smith +  pc - the preconditioner context
10069dcbbd2bSBarry Smith -  form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
10079dcbbd2bSBarry Smith    PC_MG_FULL, PC_MG_KASKADE
10084b9ad928SBarry Smith 
10094b9ad928SBarry Smith    Options Database Key:
10104b9ad928SBarry Smith .  -pc_mg_type <form> - Sets <form>, one of multiplicative,
10114b9ad928SBarry Smith    additive, full, kaskade
10124b9ad928SBarry Smith 
10134b9ad928SBarry Smith    Level: advanced
10144b9ad928SBarry Smith 
101597177400SBarry Smith .seealso: PCMGSetLevels()
10164b9ad928SBarry Smith @*/
10177087cfbeSBarry Smith PetscErrorCode  PCMGSetType(PC pc,PCMGType form)
10184b9ad928SBarry Smith {
1019f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
10204b9ad928SBarry Smith 
10214b9ad928SBarry Smith   PetscFunctionBegin;
10220700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1023c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,form,2);
102431567311SBarry Smith   mg->am = form;
10259dcbbd2bSBarry Smith   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1026c60c7ad4SBarry Smith   else pc->ops->applyrichardson = NULL;
1027c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1028c60c7ad4SBarry Smith }
1029c60c7ad4SBarry Smith 
1030c60c7ad4SBarry Smith /*@
1031c60c7ad4SBarry Smith    PCMGGetType - Determines the form of multigrid to use:
1032c60c7ad4SBarry Smith    multiplicative, additive, full, or the Kaskade algorithm.
1033c60c7ad4SBarry Smith 
1034c60c7ad4SBarry Smith    Logically Collective on PC
1035c60c7ad4SBarry Smith 
1036c60c7ad4SBarry Smith    Input Parameter:
1037c60c7ad4SBarry Smith .  pc - the preconditioner context
1038c60c7ad4SBarry Smith 
1039c60c7ad4SBarry Smith    Output Parameter:
1040c60c7ad4SBarry Smith .  type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
1041c60c7ad4SBarry Smith 
1042c60c7ad4SBarry Smith 
1043c60c7ad4SBarry Smith    Level: advanced
1044c60c7ad4SBarry Smith 
1045c60c7ad4SBarry Smith .seealso: PCMGSetLevels()
1046c60c7ad4SBarry Smith @*/
1047c60c7ad4SBarry Smith PetscErrorCode  PCMGGetType(PC pc,PCMGType *type)
1048c60c7ad4SBarry Smith {
1049c60c7ad4SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1050c60c7ad4SBarry Smith 
1051c60c7ad4SBarry Smith   PetscFunctionBegin;
1052c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1053c60c7ad4SBarry Smith   *type = mg->am;
10544b9ad928SBarry Smith   PetscFunctionReturn(0);
10554b9ad928SBarry Smith }
10564b9ad928SBarry Smith 
10574b9ad928SBarry Smith /*@
10580d353602SBarry Smith    PCMGSetCycleType - Sets the type cycles to use.  Use PCMGSetCycleTypeOnLevel() for more
10594b9ad928SBarry Smith    complicated cycling.
10604b9ad928SBarry Smith 
1061ad4df100SBarry Smith    Logically Collective on PC
10624b9ad928SBarry Smith 
10634b9ad928SBarry Smith    Input Parameters:
1064c2be2410SBarry Smith +  pc - the multigrid context
1065c1cbb1deSBarry Smith -  n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
10664b9ad928SBarry Smith 
10674b9ad928SBarry Smith    Options Database Key:
1068c1cbb1deSBarry Smith .  -pc_mg_cycle_type <v,w> - provide the cycle desired
10694b9ad928SBarry Smith 
10704b9ad928SBarry Smith    Level: advanced
10714b9ad928SBarry Smith 
10720d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel()
10734b9ad928SBarry Smith @*/
10747087cfbeSBarry Smith PetscErrorCode  PCMGSetCycleType(PC pc,PCMGCycleType n)
10754b9ad928SBarry Smith {
1076f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
1077f3fbd535SBarry Smith   PC_MG_Levels **mglevels = mg->levels;
107879416396SBarry Smith   PetscInt     i,levels;
10794b9ad928SBarry Smith 
10804b9ad928SBarry Smith   PetscFunctionBegin;
10810700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
108269fbec6eSBarry Smith   PetscValidLogicalCollectiveEnum(pc,n,2);
10831a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1084f3fbd535SBarry Smith   levels = mglevels[0]->levels;
10852fa5cd67SKarl Rupp   for (i=0; i<levels; i++) mglevels[i]->cycles = n;
10864b9ad928SBarry Smith   PetscFunctionReturn(0);
10874b9ad928SBarry Smith }
10884b9ad928SBarry Smith 
10898cc2d5dfSBarry Smith /*@
10908cc2d5dfSBarry Smith    PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
10918cc2d5dfSBarry Smith          of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
10928cc2d5dfSBarry Smith 
1093ad4df100SBarry Smith    Logically Collective on PC
10948cc2d5dfSBarry Smith 
10958cc2d5dfSBarry Smith    Input Parameters:
10968cc2d5dfSBarry Smith +  pc - the multigrid context
10978cc2d5dfSBarry Smith -  n - number of cycles (default is 1)
10988cc2d5dfSBarry Smith 
10998cc2d5dfSBarry Smith    Options Database Key:
1100e1bc860dSBarry Smith .  -pc_mg_multiplicative_cycles n
11018cc2d5dfSBarry Smith 
11028cc2d5dfSBarry Smith    Level: advanced
11038cc2d5dfSBarry Smith 
110495452b02SPatrick Sanan    Notes:
110595452b02SPatrick Sanan     This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
11068cc2d5dfSBarry Smith 
11078cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
11088cc2d5dfSBarry Smith @*/
11097087cfbeSBarry Smith PetscErrorCode  PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
11108cc2d5dfSBarry Smith {
1111f3fbd535SBarry Smith   PC_MG        *mg        = (PC_MG*)pc->data;
11128cc2d5dfSBarry Smith 
11138cc2d5dfSBarry Smith   PetscFunctionBegin;
11140700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1115c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
11162a21e185SBarry Smith   mg->cyclesperpcapply = n;
11178cc2d5dfSBarry Smith   PetscFunctionReturn(0);
11188cc2d5dfSBarry Smith }
11198cc2d5dfSBarry Smith 
11202134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1121fb15c04eSBarry Smith {
1122fb15c04eSBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
1123fb15c04eSBarry Smith 
1124fb15c04eSBarry Smith   PetscFunctionBegin;
11252134b1e4SBarry Smith   mg->galerkin = use;
1126fb15c04eSBarry Smith   PetscFunctionReturn(0);
1127fb15c04eSBarry Smith }
1128fb15c04eSBarry Smith 
1129c2be2410SBarry Smith /*@
113097177400SBarry Smith    PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
113182b87d87SMatthew G. Knepley       finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1132c2be2410SBarry Smith 
1133ad4df100SBarry Smith    Logically Collective on PC
1134c2be2410SBarry Smith 
1135c2be2410SBarry Smith    Input Parameters:
1136c91913d3SJed Brown +  pc - the multigrid context
11372134b1e4SBarry Smith -  use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1138c2be2410SBarry Smith 
1139c2be2410SBarry Smith    Options Database Key:
11402134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none>
1141c2be2410SBarry Smith 
1142c2be2410SBarry Smith    Level: intermediate
1143c2be2410SBarry Smith 
114495452b02SPatrick Sanan    Notes:
114595452b02SPatrick Sanan     Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
11461c1aac46SBarry Smith      use the PCMG construction of the coarser grids.
11471c1aac46SBarry Smith 
11482134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType
11493fc8bf9cSBarry Smith 
1150c2be2410SBarry Smith @*/
11512134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1152c2be2410SBarry Smith {
1153fb15c04eSBarry Smith   PetscErrorCode ierr;
1154c2be2410SBarry Smith 
1155c2be2410SBarry Smith   PetscFunctionBegin;
11560700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11572134b1e4SBarry Smith   ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr);
1158c2be2410SBarry Smith   PetscFunctionReturn(0);
1159c2be2410SBarry Smith }
1160c2be2410SBarry Smith 
11613fc8bf9cSBarry Smith /*@
11623fc8bf9cSBarry Smith    PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
116382b87d87SMatthew G. Knepley       A_i-1 = r_i * A_i * p_i
11643fc8bf9cSBarry Smith 
11653fc8bf9cSBarry Smith    Not Collective
11663fc8bf9cSBarry Smith 
11673fc8bf9cSBarry Smith    Input Parameter:
11683fc8bf9cSBarry Smith .  pc - the multigrid context
11693fc8bf9cSBarry Smith 
11703fc8bf9cSBarry Smith    Output Parameter:
11712134b1e4SBarry 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
11723fc8bf9cSBarry Smith 
11733fc8bf9cSBarry Smith    Level: intermediate
11743fc8bf9cSBarry Smith 
11752134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType
11763fc8bf9cSBarry Smith 
11773fc8bf9cSBarry Smith @*/
11782134b1e4SBarry Smith PetscErrorCode  PCMGGetGalerkin(PC pc,PCMGGalerkinType  *galerkin)
11793fc8bf9cSBarry Smith {
1180f3fbd535SBarry Smith   PC_MG *mg = (PC_MG*)pc->data;
11813fc8bf9cSBarry Smith 
11823fc8bf9cSBarry Smith   PetscFunctionBegin;
11830700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
11842134b1e4SBarry Smith   *galerkin = mg->galerkin;
11853fc8bf9cSBarry Smith   PetscFunctionReturn(0);
11863fc8bf9cSBarry Smith }
11873fc8bf9cSBarry Smith 
1188f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1189f3b08a26SMatthew G. Knepley {
1190f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
1191f3b08a26SMatthew G. Knepley 
1192f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1193f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = adapt;
1194f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1195f3b08a26SMatthew G. Knepley }
1196f3b08a26SMatthew G. Knepley 
1197f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1198f3b08a26SMatthew G. Knepley {
1199f3b08a26SMatthew G. Knepley   PC_MG *mg = (PC_MG *) pc->data;
1200f3b08a26SMatthew G. Knepley 
1201f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1202f3b08a26SMatthew G. Knepley   *adapt = mg->adaptInterpolation;
1203f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1204f3b08a26SMatthew G. Knepley }
1205f3b08a26SMatthew G. Knepley 
1206f3b08a26SMatthew G. Knepley /*@
1207f3b08a26SMatthew G. Knepley   PCMGSetAdaptInterpolation - Adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated.
1208f3b08a26SMatthew G. Knepley 
1209f3b08a26SMatthew G. Knepley   Logically Collective on PC
1210f3b08a26SMatthew G. Knepley 
1211f3b08a26SMatthew G. Knepley   Input Parameters:
1212f3b08a26SMatthew G. Knepley + pc    - the multigrid context
1213f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator
1214f3b08a26SMatthew G. Knepley 
1215f3b08a26SMatthew G. Knepley   Options Database Keys:
1216f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp                     - Turn on adaptation
1217f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int>             - The number of modes to use, should be divisible by dimension
1218f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector
1219f3b08a26SMatthew G. Knepley 
1220f3b08a26SMatthew G. Knepley   Level: intermediate
1221f3b08a26SMatthew G. Knepley 
1222f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin
1223f3b08a26SMatthew G. Knepley .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1224f3b08a26SMatthew G. Knepley @*/
1225f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1226f3b08a26SMatthew G. Knepley {
1227f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1228f3b08a26SMatthew G. Knepley 
1229f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1230f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1231f3b08a26SMatthew G. Knepley   ierr = PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));CHKERRQ(ierr);
1232f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1233f3b08a26SMatthew G. Knepley }
1234f3b08a26SMatthew G. Knepley 
1235f3b08a26SMatthew G. Knepley /*@
1236f3b08a26SMatthew G. Knepley   PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated.
1237f3b08a26SMatthew G. Knepley 
1238f3b08a26SMatthew G. Knepley   Not collective
1239f3b08a26SMatthew G. Knepley 
1240f3b08a26SMatthew G. Knepley   Input Parameter:
1241f3b08a26SMatthew G. Knepley . pc    - the multigrid context
1242f3b08a26SMatthew G. Knepley 
1243f3b08a26SMatthew G. Knepley   Output Parameter:
1244f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator
1245f3b08a26SMatthew G. Knepley 
1246f3b08a26SMatthew G. Knepley   Level: intermediate
1247f3b08a26SMatthew G. Knepley 
1248f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin
1249f3b08a26SMatthew G. Knepley .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1250f3b08a26SMatthew G. Knepley @*/
1251f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1252f3b08a26SMatthew G. Knepley {
1253f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1254f3b08a26SMatthew G. Knepley 
1255f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1256f3b08a26SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1257f3b08a26SMatthew G. Knepley   PetscValidPointer(adapt, 2);
1258f3b08a26SMatthew G. Knepley   ierr = PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));CHKERRQ(ierr);
1259f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1260f3b08a26SMatthew G. Knepley }
1261f3b08a26SMatthew G. Knepley 
12624b9ad928SBarry Smith /*@
126306792cafSBarry Smith    PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1264710315b6SLawrence Mitchell    on all levels.  Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1265710315b6SLawrence Mitchell    pre- and post-smoothing steps.
126606792cafSBarry Smith 
126706792cafSBarry Smith    Logically Collective on PC
126806792cafSBarry Smith 
126906792cafSBarry Smith    Input Parameters:
127006792cafSBarry Smith +  mg - the multigrid context
127106792cafSBarry Smith -  n - the number of smoothing steps
127206792cafSBarry Smith 
127306792cafSBarry Smith    Options Database Key:
1274a2b725a8SWilliam Gropp .  -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
127506792cafSBarry Smith 
127606792cafSBarry Smith    Level: advanced
127706792cafSBarry Smith 
127895452b02SPatrick Sanan    Notes:
127995452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
128006792cafSBarry Smith     there is no separate smooth up on the coarsest grid.
128106792cafSBarry Smith 
1282710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp()
128306792cafSBarry Smith @*/
128406792cafSBarry Smith PetscErrorCode  PCMGSetNumberSmooth(PC pc,PetscInt n)
128506792cafSBarry Smith {
128606792cafSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
128706792cafSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
128806792cafSBarry Smith   PetscErrorCode ierr;
128906792cafSBarry Smith   PetscInt       i,levels;
129006792cafSBarry Smith 
129106792cafSBarry Smith   PetscFunctionBegin;
129206792cafSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
129306792cafSBarry Smith   PetscValidLogicalCollectiveInt(pc,n,2);
12941a97d7b6SStefano Zampini   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
129506792cafSBarry Smith   levels = mglevels[0]->levels;
129606792cafSBarry Smith 
129706792cafSBarry Smith   for (i=1; i<levels; i++) {
129806792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
129906792cafSBarry Smith     ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr);
130006792cafSBarry Smith     mg->default_smoothu = n;
130106792cafSBarry Smith     mg->default_smoothd = n;
130206792cafSBarry Smith   }
130306792cafSBarry Smith   PetscFunctionReturn(0);
130406792cafSBarry Smith }
130506792cafSBarry Smith 
1306f442ab6aSBarry Smith /*@
13078e5aa403SBarry Smith    PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels
1308710315b6SLawrence Mitchell        and adds the suffix _up to the options name
1309f442ab6aSBarry Smith 
1310f442ab6aSBarry Smith    Logically Collective on PC
1311f442ab6aSBarry Smith 
1312f442ab6aSBarry Smith    Input Parameters:
1313f442ab6aSBarry Smith .  pc - the preconditioner context
1314f442ab6aSBarry Smith 
1315f442ab6aSBarry Smith    Options Database Key:
1316f442ab6aSBarry Smith .  -pc_mg_distinct_smoothup
1317f442ab6aSBarry Smith 
1318f442ab6aSBarry Smith    Level: advanced
1319f442ab6aSBarry Smith 
132095452b02SPatrick Sanan    Notes:
132195452b02SPatrick Sanan     this does not set a value on the coarsest grid, since we assume that
1322f442ab6aSBarry Smith     there is no separate smooth up on the coarsest grid.
1323f442ab6aSBarry Smith 
1324710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth()
1325f442ab6aSBarry Smith @*/
1326f442ab6aSBarry Smith PetscErrorCode  PCMGSetDistinctSmoothUp(PC pc)
1327f442ab6aSBarry Smith {
1328f442ab6aSBarry Smith   PC_MG          *mg        = (PC_MG*)pc->data;
1329f442ab6aSBarry Smith   PC_MG_Levels   **mglevels = mg->levels;
1330f442ab6aSBarry Smith   PetscErrorCode ierr;
1331f442ab6aSBarry Smith   PetscInt       i,levels;
1332f442ab6aSBarry Smith   KSP            subksp;
1333f442ab6aSBarry Smith 
1334f442ab6aSBarry Smith   PetscFunctionBegin;
1335f442ab6aSBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
1336f442ab6aSBarry Smith   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1337f442ab6aSBarry Smith   levels = mglevels[0]->levels;
1338f442ab6aSBarry Smith 
1339f442ab6aSBarry Smith   for (i=1; i<levels; i++) {
1340710315b6SLawrence Mitchell     const char *prefix = NULL;
1341f442ab6aSBarry Smith     /* make sure smoother up and down are different */
1342f442ab6aSBarry Smith     ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr);
1343710315b6SLawrence Mitchell     ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr);
1344710315b6SLawrence Mitchell     ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr);
1345f442ab6aSBarry Smith     ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr);
1346f442ab6aSBarry Smith   }
1347f442ab6aSBarry Smith   PetscFunctionReturn(0);
1348f442ab6aSBarry Smith }
1349f442ab6aSBarry Smith 
135007a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
135107a4832bSFande Kong PetscErrorCode  PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[])
135207a4832bSFande Kong {
135307a4832bSFande Kong   PC_MG          *mg        = (PC_MG*)pc->data;
135407a4832bSFande Kong   PC_MG_Levels   **mglevels = mg->levels;
135507a4832bSFande Kong   Mat            *mat;
135607a4832bSFande Kong   PetscInt       l;
135707a4832bSFande Kong   PetscErrorCode ierr;
135807a4832bSFande Kong 
135907a4832bSFande Kong   PetscFunctionBegin;
136007a4832bSFande Kong   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
136107a4832bSFande Kong   ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr);
136207a4832bSFande Kong   for (l=1; l< mg->nlevels; l++) {
136307a4832bSFande Kong     mat[l-1] = mglevels[l]->interpolate;
136407a4832bSFande Kong     ierr = PetscObjectReference((PetscObject)mat[l-1]);CHKERRQ(ierr);
136507a4832bSFande Kong   }
136607a4832bSFande Kong   *num_levels = mg->nlevels;
136707a4832bSFande Kong   *interpolations = mat;
136807a4832bSFande Kong   PetscFunctionReturn(0);
136907a4832bSFande Kong }
137007a4832bSFande Kong 
137107a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
137207a4832bSFande Kong PetscErrorCode  PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[])
137307a4832bSFande Kong {
137407a4832bSFande Kong   PC_MG          *mg        = (PC_MG*)pc->data;
137507a4832bSFande Kong   PC_MG_Levels   **mglevels = mg->levels;
137607a4832bSFande Kong   PetscInt       l;
137707a4832bSFande Kong   Mat            *mat;
137807a4832bSFande Kong   PetscErrorCode ierr;
137907a4832bSFande Kong 
138007a4832bSFande Kong   PetscFunctionBegin;
138107a4832bSFande Kong   if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
138207a4832bSFande Kong   ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr);
138307a4832bSFande Kong   for (l=0; l<mg->nlevels-1; l++) {
138407a4832bSFande Kong     ierr = KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));CHKERRQ(ierr);
138507a4832bSFande Kong     ierr = PetscObjectReference((PetscObject)mat[l]);CHKERRQ(ierr);
138607a4832bSFande Kong   }
138707a4832bSFande Kong   *num_levels = mg->nlevels;
138807a4832bSFande Kong   *coarseOperators = mat;
138907a4832bSFande Kong   PetscFunctionReturn(0);
139007a4832bSFande Kong }
139107a4832bSFande Kong 
1392f3b08a26SMatthew G. Knepley /*@C
1393f3b08a26SMatthew G. Knepley   PCMGRegisterCoarseSpaceConstructor -  Adds a method to the PCMG package for coarse space construction.
1394f3b08a26SMatthew G. Knepley 
1395f3b08a26SMatthew G. Knepley   Not collective
1396f3b08a26SMatthew G. Knepley 
1397f3b08a26SMatthew G. Knepley   Input Parameters:
1398f3b08a26SMatthew G. Knepley + name     - name of the constructor
1399f3b08a26SMatthew G. Knepley - function - constructor routine
1400f3b08a26SMatthew G. Knepley 
1401f3b08a26SMatthew G. Knepley   Notes:
1402f3b08a26SMatthew G. Knepley   Calling sequence for the routine:
1403f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1404f3b08a26SMatthew G. Knepley $   pc        - The PC object
1405f3b08a26SMatthew G. Knepley $   l         - The multigrid level, 0 is the coarse level
1406f3b08a26SMatthew G. Knepley $   dm        - The DM for this level
1407f3b08a26SMatthew G. Knepley $   smooth    - The level smoother
1408f3b08a26SMatthew G. Knepley $   Nc        - The size of the coarse space
1409f3b08a26SMatthew G. Knepley $   initGuess - Basis for an initial guess for the space
1410f3b08a26SMatthew G. Knepley $   coarseSp  - A basis for the computed coarse space
1411f3b08a26SMatthew G. Knepley 
1412f3b08a26SMatthew G. Knepley   Level: advanced
1413f3b08a26SMatthew G. Knepley 
1414f3b08a26SMatthew G. Knepley .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister()
1415f3b08a26SMatthew G. Knepley @*/
1416f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1417f3b08a26SMatthew G. Knepley {
1418f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1419f3b08a26SMatthew G. Knepley 
1420f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1421f3b08a26SMatthew G. Knepley   ierr = PCInitializePackage();CHKERRQ(ierr);
1422f3b08a26SMatthew G. Knepley   ierr = PetscFunctionListAdd(&PCMGCoarseList,name,function);CHKERRQ(ierr);
1423f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1424f3b08a26SMatthew G. Knepley }
1425f3b08a26SMatthew G. Knepley 
1426f3b08a26SMatthew G. Knepley /*@C
1427f3b08a26SMatthew G. Knepley   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.
1428f3b08a26SMatthew G. Knepley 
1429f3b08a26SMatthew G. Knepley   Not collective
1430f3b08a26SMatthew G. Knepley 
1431f3b08a26SMatthew G. Knepley   Input Parameter:
1432f3b08a26SMatthew G. Knepley . name     - name of the constructor
1433f3b08a26SMatthew G. Knepley 
1434f3b08a26SMatthew G. Knepley   Output Parameter:
1435f3b08a26SMatthew G. Knepley . function - constructor routine
1436f3b08a26SMatthew G. Knepley 
1437f3b08a26SMatthew G. Knepley   Notes:
1438f3b08a26SMatthew G. Knepley   Calling sequence for the routine:
1439f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1440f3b08a26SMatthew G. Knepley $   pc        - The PC object
1441f3b08a26SMatthew G. Knepley $   l         - The multigrid level, 0 is the coarse level
1442f3b08a26SMatthew G. Knepley $   dm        - The DM for this level
1443f3b08a26SMatthew G. Knepley $   smooth    - The level smoother
1444f3b08a26SMatthew G. Knepley $   Nc        - The size of the coarse space
1445f3b08a26SMatthew G. Knepley $   initGuess - Basis for an initial guess for the space
1446f3b08a26SMatthew G. Knepley $   coarseSp  - A basis for the computed coarse space
1447f3b08a26SMatthew G. Knepley 
1448f3b08a26SMatthew G. Knepley   Level: advanced
1449f3b08a26SMatthew G. Knepley 
1450f3b08a26SMatthew G. Knepley .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister()
1451f3b08a26SMatthew G. Knepley @*/
1452f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1453f3b08a26SMatthew G. Knepley {
1454f3b08a26SMatthew G. Knepley   PetscErrorCode ierr;
1455f3b08a26SMatthew G. Knepley 
1456f3b08a26SMatthew G. Knepley   PetscFunctionBegin;
1457f3b08a26SMatthew G. Knepley   ierr = PetscFunctionListFind(PCMGCoarseList,name,function);CHKERRQ(ierr);
1458f3b08a26SMatthew G. Knepley   PetscFunctionReturn(0);
1459f3b08a26SMatthew G. Knepley }
1460f3b08a26SMatthew G. Knepley 
14614b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/
14624b9ad928SBarry Smith 
14633b09bd56SBarry Smith /*MC
1464ccb205f8SBarry Smith    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
14653b09bd56SBarry Smith     information about the coarser grid matrices and restriction/interpolation operators.
14663b09bd56SBarry Smith 
14673b09bd56SBarry Smith    Options Database Keys:
14683b09bd56SBarry Smith +  -pc_mg_levels <nlevels> - number of levels including finest
1469391689abSStefano Zampini .  -pc_mg_cycle_type <v,w> - provide the cycle desired
14708c1c2452SJed Brown .  -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
14713b09bd56SBarry Smith .  -pc_mg_log - log information about time spent on each level of the solver
1472710315b6SLawrence Mitchell .  -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
14732134b1e4SBarry Smith .  -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
14748cf6037eSBarry Smith .  -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
14758cf6037eSBarry Smith .  -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1476e3c5b3baSBarry Smith                         to the Socket viewer for reading from MATLAB.
14778cf6037eSBarry Smith -  -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
14788cf6037eSBarry Smith                         to the binary output file called binaryoutput
14793b09bd56SBarry Smith 
148095452b02SPatrick Sanan    Notes:
14810b3c753eSRichard 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
14823b09bd56SBarry Smith 
14838cf6037eSBarry Smith        When run with a single level the smoother options are used on that level NOT the coarse grid solver options
14848cf6037eSBarry Smith 
148523067569SBarry Smith        When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
148623067569SBarry Smith        is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
148723067569SBarry Smith        (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
148823067569SBarry Smith        residual is computed at the end of each cycle.
148923067569SBarry Smith 
14903b09bd56SBarry Smith    Level: intermediate
14913b09bd56SBarry Smith 
14928cf6037eSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1493710315b6SLawrence Mitchell            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1494710315b6SLawrence Mitchell            PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
149597177400SBarry Smith            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
14960d353602SBarry Smith            PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
14973b09bd56SBarry Smith M*/
14983b09bd56SBarry Smith 
14998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
15004b9ad928SBarry Smith {
1501f3fbd535SBarry Smith   PC_MG          *mg;
1502f3fbd535SBarry Smith   PetscErrorCode ierr;
1503f3fbd535SBarry Smith 
15044b9ad928SBarry Smith   PetscFunctionBegin;
1505b00a9115SJed Brown   ierr         = PetscNewLog(pc,&mg);CHKERRQ(ierr);
1506f3fbd535SBarry Smith   pc->data     = (void*)mg;
1507f3fbd535SBarry Smith   mg->nlevels  = -1;
150810eca3edSLisandro Dalcin   mg->am       = PC_MG_MULTIPLICATIVE;
15092134b1e4SBarry Smith   mg->galerkin = PC_MG_GALERKIN_NONE;
1510f3b08a26SMatthew G. Knepley   mg->adaptInterpolation = PETSC_FALSE;
1511f3b08a26SMatthew G. Knepley   mg->Nc                 = -1;
1512f3b08a26SMatthew G. Knepley   mg->eigenvalue         = -1;
1513f3fbd535SBarry Smith 
151437a44384SMark Adams   pc->useAmat = PETSC_TRUE;
151537a44384SMark Adams 
15164b9ad928SBarry Smith   pc->ops->apply          = PCApply_MG;
15174b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_MG;
1518a06653b4SBarry Smith   pc->ops->reset          = PCReset_MG;
15194b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_MG;
15204b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_MG;
15214b9ad928SBarry Smith   pc->ops->view           = PCView_MG;
1522fb15c04eSBarry Smith 
1523f3b08a26SMatthew G. Knepley   ierr = PetscObjectComposedDataRegister(&mg->eigenvalue);CHKERRQ(ierr);
1524fb15c04eSBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr);
152597d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr);
152697d33e41SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr);
1527fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);CHKERRQ(ierr);
1528fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);CHKERRQ(ierr);
1529f3b08a26SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);CHKERRQ(ierr);
1530f3b08a26SMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);CHKERRQ(ierr);
15314b9ad928SBarry Smith   PetscFunctionReturn(0);
15324b9ad928SBarry Smith }
1533