Actual source code: mg.c
petsc-3.4.5 2014-06-29
2: /*
3: Defines the multigrid preconditioner interface.
4: */
5: #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscksp.h" I*/
6: #include <petscdm.h>
10: PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason)
11: {
12: PC_MG *mg = (PC_MG*)pc->data;
13: PC_MG_Levels *mgc,*mglevels = *mglevelsin;
15: PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
18: if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
19: KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x); /* pre-smooth */
20: if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
21: if (mglevels->level) { /* not the coarsest grid */
22: if (mglevels->eventresidual) {PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);}
23: (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);
24: if (mglevels->eventresidual) {PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);}
26: /* if on finest level and have convergence criteria set */
27: if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
28: PetscReal rnorm;
29: VecNorm(mglevels->r,NORM_2,&rnorm);
30: if (rnorm <= mg->ttol) {
31: if (rnorm < mg->abstol) {
32: *reason = PCRICHARDSON_CONVERGED_ATOL;
33: PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);
34: } else {
35: *reason = PCRICHARDSON_CONVERGED_RTOL;
36: PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);
37: }
38: return(0);
39: }
40: }
42: mgc = *(mglevelsin - 1);
43: if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
44: MatRestrict(mglevels->restrct,mglevels->r,mgc->b);
45: if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
46: VecSet(mgc->x,0.0);
47: while (cycles--) {
48: PCMGMCycle_Private(pc,mglevelsin-1,reason);
49: }
50: if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
51: MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);
52: if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
53: if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
54: KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x); /* post smooth */
55: if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
56: }
57: return(0);
58: }
62: 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)
63: {
64: PC_MG *mg = (PC_MG*)pc->data;
65: PC_MG_Levels **mglevels = mg->levels;
67: PetscInt levels = mglevels[0]->levels,i;
70: mglevels[levels-1]->b = b;
71: mglevels[levels-1]->x = x;
73: mg->rtol = rtol;
74: mg->abstol = abstol;
75: mg->dtol = dtol;
76: if (rtol) {
77: /* compute initial residual norm for relative convergence test */
78: PetscReal rnorm;
79: if (zeroguess) {
80: VecNorm(b,NORM_2,&rnorm);
81: } else {
82: (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);
83: VecNorm(w,NORM_2,&rnorm);
84: }
85: mg->ttol = PetscMax(rtol*rnorm,abstol);
86: } else if (abstol) mg->ttol = abstol;
87: else mg->ttol = 0.0;
89: /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
90: stop prematurely due to small residual */
91: for (i=1; i<levels; i++) {
92: KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
93: if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
94: KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
95: }
96: }
98: *reason = (PCRichardsonConvergedReason)0;
99: for (i=0; i<its; i++) {
100: PCMGMCycle_Private(pc,mglevels+levels-1,reason);
101: if (*reason) break;
102: }
103: if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
104: *outits = i;
105: return(0);
106: }
110: PetscErrorCode PCReset_MG(PC pc)
111: {
112: PC_MG *mg = (PC_MG*)pc->data;
113: PC_MG_Levels **mglevels = mg->levels;
115: PetscInt i,n;
118: if (mglevels) {
119: n = mglevels[0]->levels;
120: for (i=0; i<n-1; i++) {
121: VecDestroy(&mglevels[i+1]->r);
122: VecDestroy(&mglevels[i]->b);
123: VecDestroy(&mglevels[i]->x);
124: MatDestroy(&mglevels[i+1]->restrct);
125: MatDestroy(&mglevels[i+1]->interpolate);
126: VecDestroy(&mglevels[i+1]->rscale);
127: }
129: for (i=0; i<n; i++) {
130: MatDestroy(&mglevels[i]->A);
131: if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
132: KSPReset(mglevels[i]->smoothd);
133: }
134: KSPReset(mglevels[i]->smoothu);
135: }
136: }
137: return(0);
138: }
142: /*@C
143: PCMGSetLevels - Sets the number of levels to use with MG.
144: Must be called before any other MG routine.
146: Logically Collective on PC
148: Input Parameters:
149: + pc - the preconditioner context
150: . levels - the number of levels
151: - comms - optional communicators for each level; this is to allow solving the coarser problems
152: on smaller sets of processors. Use NULL_OBJECT for default in Fortran
154: Level: intermediate
156: Notes:
157: If the number of levels is one then the multigrid uses the -mg_levels prefix
158: for setting the level options rather than the -mg_coarse prefix.
160: .keywords: MG, set, levels, multigrid
162: .seealso: PCMGSetType(), PCMGGetLevels()
163: @*/
164: PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
165: {
167: PC_MG *mg = (PC_MG*)pc->data;
168: MPI_Comm comm;
169: PC_MG_Levels **mglevels = mg->levels;
170: PetscInt i;
171: PetscMPIInt size;
172: const char *prefix;
173: PC ipc;
174: PetscInt n;
179: PetscObjectGetComm((PetscObject)pc,&comm);
180: if (mg->nlevels == levels) return(0);
181: if (mglevels) {
182: /* changing the number of levels so free up the previous stuff */
183: PCReset_MG(pc);
184: n = mglevels[0]->levels;
185: for (i=0; i<n; i++) {
186: if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
187: KSPDestroy(&mglevels[i]->smoothd);
188: }
189: KSPDestroy(&mglevels[i]->smoothu);
190: PetscFree(mglevels[i]);
191: }
192: PetscFree(mg->levels);
193: }
195: mg->nlevels = levels;
197: PetscMalloc(levels*sizeof(PC_MG*),&mglevels);
198: PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));
200: PCGetOptionsPrefix(pc,&prefix);
202: mg->stageApply = 0;
203: for (i=0; i<levels; i++) {
204: PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);
206: mglevels[i]->level = i;
207: mglevels[i]->levels = levels;
208: mglevels[i]->cycles = PC_MG_CYCLE_V;
209: mg->default_smoothu = 2;
210: mg->default_smoothd = 2;
211: mglevels[i]->eventsmoothsetup = 0;
212: mglevels[i]->eventsmoothsolve = 0;
213: mglevels[i]->eventresidual = 0;
214: mglevels[i]->eventinterprestrict = 0;
216: if (comms) comm = comms[i];
217: KSPCreate(comm,&mglevels[i]->smoothd);
218: KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);
219: KSPSetConvergenceTest(mglevels[i]->smoothd,KSPSkipConverged,NULL,NULL);
220: KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);
221: KSPGetPC(mglevels[i]->smoothd,&ipc);
222: PCSetType(ipc,PCSOR);
223: PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);
224: KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, i ? mg->default_smoothd : 1);
225: KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);
227: /* do special stuff for coarse grid */
228: if (!i && levels > 1) {
229: KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");
231: /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
232: KSPSetType(mglevels[0]->smoothd,KSPPREONLY);
233: KSPGetPC(mglevels[0]->smoothd,&ipc);
234: MPI_Comm_size(comm,&size);
235: if (size > 1) {
236: KSP innerksp;
237: PC innerpc;
238: PCSetType(ipc,PCREDUNDANT);
239: PCRedundantGetKSP(ipc,&innerksp);
240: KSPGetPC(innerksp,&innerpc);
241: PCFactorSetShiftType(innerpc,MAT_SHIFT_INBLOCKS);
242: } else {
243: PCSetType(ipc,PCLU);
244: PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);
245: }
246: } else {
247: char tprefix[128];
248: sprintf(tprefix,"mg_levels_%d_",(int)i);
249: KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);
250: }
251: PetscLogObjectParent(pc,mglevels[i]->smoothd);
253: mglevels[i]->smoothu = mglevels[i]->smoothd;
254: mg->rtol = 0.0;
255: mg->abstol = 0.0;
256: mg->dtol = 0.0;
257: mg->ttol = 0.0;
258: mg->cyclesperpcapply = 1;
259: }
260: mg->am = PC_MG_MULTIPLICATIVE;
261: mg->levels = mglevels;
262: pc->ops->applyrichardson = PCApplyRichardson_MG;
263: return(0);
264: }
269: PetscErrorCode PCDestroy_MG(PC pc)
270: {
272: PC_MG *mg = (PC_MG*)pc->data;
273: PC_MG_Levels **mglevels = mg->levels;
274: PetscInt i,n;
277: PCReset_MG(pc);
278: if (mglevels) {
279: n = mglevels[0]->levels;
280: for (i=0; i<n; i++) {
281: if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
282: KSPDestroy(&mglevels[i]->smoothd);
283: }
284: KSPDestroy(&mglevels[i]->smoothu);
285: PetscFree(mglevels[i]);
286: }
287: PetscFree(mg->levels);
288: }
289: PetscFree(pc->data);
290: return(0);
291: }
295: extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**);
296: extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**);
297: extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**);
299: /*
300: PCApply_MG - Runs either an additive, multiplicative, Kaskadic
301: or full cycle of multigrid.
303: Note:
304: A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
305: */
308: static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
309: {
310: PC_MG *mg = (PC_MG*)pc->data;
311: PC_MG_Levels **mglevels = mg->levels;
313: PetscInt levels = mglevels[0]->levels,i;
316: if (mg->stageApply) {PetscLogStagePush(mg->stageApply);}
317: /* When the DM is supplying the matrix then it will not exist until here */
318: for (i=0; i<levels; i++) {
319: if (!mglevels[i]->A) {
320: KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL,NULL);
321: PetscObjectReference((PetscObject)mglevels[i]->A);
322: }
323: }
325: mglevels[levels-1]->b = b;
326: mglevels[levels-1]->x = x;
327: if (mg->am == PC_MG_MULTIPLICATIVE) {
328: VecSet(x,0.0);
329: for (i=0; i<mg->cyclesperpcapply; i++) {
330: PCMGMCycle_Private(pc,mglevels+levels-1,NULL);
331: }
332: } else if (mg->am == PC_MG_ADDITIVE) {
333: PCMGACycle_Private(pc,mglevels);
334: } else if (mg->am == PC_MG_KASKADE) {
335: PCMGKCycle_Private(pc,mglevels);
336: } else {
337: PCMGFCycle_Private(pc,mglevels);
338: }
339: if (mg->stageApply) {PetscLogStagePop();}
340: return(0);
341: }
346: PetscErrorCode PCSetFromOptions_MG(PC pc)
347: {
349: PetscInt m,levels = 1,cycles;
350: PetscBool flg,set;
351: PC_MG *mg = (PC_MG*)pc->data;
352: PC_MG_Levels **mglevels = mg->levels;
353: PCMGType mgtype;
354: PCMGCycleType mgctype;
357: PetscOptionsHead("Multigrid options");
358: if (!mg->levels) {
359: PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
360: if (!flg && pc->dm) {
361: DMGetRefineLevel(pc->dm,&levels);
362: levels++;
363: mg->usedmfornumberoflevels = PETSC_TRUE;
364: }
365: PCMGSetLevels(pc,levels,NULL);
366: }
367: mglevels = mg->levels;
369: mgctype = (PCMGCycleType) mglevels[0]->cycles;
370: PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);
371: if (flg) {
372: PCMGSetCycleType(pc,mgctype);
373: }
374: flg = PETSC_FALSE;
375: PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,&set);
376: if (set) {
377: PCMGSetGalerkin(pc,flg);
378: }
379: PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);
380: if (flg) {
381: PCMGSetNumberSmoothUp(pc,m);
382: }
383: PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);
384: if (flg) {
385: PCMGSetNumberSmoothDown(pc,m);
386: }
387: mgtype = mg->am;
388: PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);
389: if (flg) {
390: PCMGSetType(pc,mgtype);
391: }
392: if (mg->am == PC_MG_MULTIPLICATIVE) {
393: PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);
394: if (flg) {
395: PCMGMultiplicativeSetCycles(pc,cycles);
396: }
397: }
398: flg = PETSC_FALSE;
399: PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);
400: if (flg) {
401: PetscInt i;
402: char eventname[128];
403: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
404: levels = mglevels[0]->levels;
405: for (i=0; i<levels; i++) {
406: sprintf(eventname,"MGSetup Level %d",(int)i);
407: PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);
408: sprintf(eventname,"MGSmooth Level %d",(int)i);
409: PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);
410: if (i) {
411: sprintf(eventname,"MGResid Level %d",(int)i);
412: PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);
413: sprintf(eventname,"MGInterp Level %d",(int)i);
414: PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);
415: }
416: }
418: #if defined(PETSC_USE_LOG)
419: {
420: const char *sname = "MG Apply";
421: PetscStageLog stageLog;
422: PetscInt st;
425: PetscLogGetStageLog(&stageLog);
426: for (st = 0; st < stageLog->numStages; ++st) {
427: PetscBool same;
429: PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);
430: if (same) mg->stageApply = st;
431: }
432: if (!mg->stageApply) {
433: PetscLogStageRegister(sname, &mg->stageApply);
434: }
435: }
436: #endif
437: }
438: PetscOptionsTail();
439: return(0);
440: }
442: const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0};
443: const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0};
445: #include <petscdraw.h>
448: PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
449: {
450: PC_MG *mg = (PC_MG*)pc->data;
451: PC_MG_Levels **mglevels = mg->levels;
453: PetscInt levels = mglevels ? mglevels[0]->levels : 0,i;
454: PetscBool iascii,isbinary,isdraw;
457: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
458: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
459: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
460: if (iascii) {
461: const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
462: PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);
463: if (mg->am == PC_MG_MULTIPLICATIVE) {
464: PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);
465: }
466: if (mg->galerkin) {
467: PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");
468: } else {
469: PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");
470: }
471: for (i=0; i<levels; i++) {
472: if (!i) {
473: PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);
474: } else {
475: PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);
476: }
477: PetscViewerASCIIPushTab(viewer);
478: KSPView(mglevels[i]->smoothd,viewer);
479: PetscViewerASCIIPopTab(viewer);
480: if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
481: PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");
482: } else if (i) {
483: PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);
484: PetscViewerASCIIPushTab(viewer);
485: KSPView(mglevels[i]->smoothu,viewer);
486: PetscViewerASCIIPopTab(viewer);
487: }
488: }
489: } else if (isbinary) {
490: for (i=levels-1; i>=0; i--) {
491: KSPView(mglevels[i]->smoothd,viewer);
492: if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
493: KSPView(mglevels[i]->smoothu,viewer);
494: }
495: }
496: } else if (isdraw) {
497: PetscDraw draw;
498: PetscReal x,w,y,bottom,th;
499: PetscViewerDrawGetDraw(viewer,0,&draw);
500: PetscDrawGetCurrentPoint(draw,&x,&y);
501: PetscDrawStringGetSize(draw,NULL,&th);
502: bottom = y - th;
503: for (i=levels-1; i>=0; i--) {
504: if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
505: PetscDrawPushCurrentPoint(draw,x,bottom);
506: KSPView(mglevels[i]->smoothd,viewer);
507: PetscDrawPopCurrentPoint(draw);
508: } else {
509: w = 0.5*PetscMin(1.0-x,x);
510: PetscDrawPushCurrentPoint(draw,x+w,bottom);
511: KSPView(mglevels[i]->smoothd,viewer);
512: PetscDrawPopCurrentPoint(draw);
513: PetscDrawPushCurrentPoint(draw,x-w,bottom);
514: KSPView(mglevels[i]->smoothu,viewer);
515: PetscDrawPopCurrentPoint(draw);
516: }
517: PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);
518: bottom -= th;
519: }
520: }
521: return(0);
522: }
524: #include <petsc-private/dmimpl.h>
525: #include <petsc-private/kspimpl.h>
527: /*
528: Calls setup for the KSP on each level
529: */
532: PetscErrorCode PCSetUp_MG(PC pc)
533: {
534: PC_MG *mg = (PC_MG*)pc->data;
535: PC_MG_Levels **mglevels = mg->levels;
537: PetscInt i,n = mglevels[0]->levels;
538: PC cpc;
539: PetscBool preonly,lu,redundant,cholesky,svd,dump = PETSC_FALSE,opsset,use_amat;
540: Mat dA,dB;
541: MatStructure uflag;
542: Vec tvec;
543: DM *dms;
544: PetscViewer viewer = 0;
547: /* FIX: Move this to PCSetFromOptions_MG? */
548: if (mg->usedmfornumberoflevels) {
549: PetscInt levels;
550: DMGetRefineLevel(pc->dm,&levels);
551: levels++;
552: if (levels > n) { /* the problem is now being solved on a finer grid */
553: PCMGSetLevels(pc,levels,NULL);
554: n = levels;
555: PCSetFromOptions(pc); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
556: mglevels = mg->levels;
557: }
558: }
559: KSPGetPC(mglevels[0]->smoothd,&cpc);
562: /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
563: /* so use those from global PC */
564: /* Is this what we always want? What if user wants to keep old one? */
565: KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);
566: if (opsset) {
567: Mat mmat;
568: KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat,NULL);
569: if (mmat == pc->pmat) opsset = PETSC_FALSE;
570: }
572: if (!opsset) {
573: PCGetUseAmat(pc,&use_amat);
574: if(use_amat){
575: PetscInfo(pc,"Using outer operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
576: KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);
577: }
578: else {
579: PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
580: KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat,pc->flag);
581: }
582: }
584: /* Skipping this for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? */
585: if (pc->dm && mg->galerkin != 2 && !pc->setupcalled) {
586: /* construct the interpolation from the DMs */
587: Mat p;
588: Vec rscale;
589: PetscMalloc(n*sizeof(DM),&dms);
590: dms[n-1] = pc->dm;
591: for (i=n-2; i>-1; i--) {
592: DMKSP kdm;
593: DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);
594: KSPSetDM(mglevels[i]->smoothd,dms[i]);
595: if (mg->galerkin) {KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);}
596: DMGetDMKSPWrite(dms[i],&kdm);
597: /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
598: * a bitwise OR of computing the matrix, RHS, and initial iterate. */
599: kdm->ops->computerhs = NULL;
600: kdm->rhsctx = NULL;
601: if (!mglevels[i+1]->interpolate) {
602: DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);
603: PCMGSetInterpolation(pc,i+1,p);
604: if (rscale) {PCMGSetRScale(pc,i+1,rscale);}
605: VecDestroy(&rscale);
606: MatDestroy(&p);
607: }
608: }
610: for (i=n-2; i>-1; i--) {
611: DMDestroy(&dms[i]);
612: }
613: PetscFree(dms);
614: }
616: if (pc->dm && !pc->setupcalled) {
617: /* finest smoother also gets DM but it is not active, independent of whether galerkin==2 */
618: KSPSetDM(mglevels[n-1]->smoothd,pc->dm);
619: KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);
620: }
622: if (mg->galerkin == 1) {
623: Mat B;
624: /* currently only handle case where mat and pmat are the same on coarser levels */
625: KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);
626: if (!pc->setupcalled) {
627: for (i=n-2; i>-1; i--) {
628: MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);
629: KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);
630: if (i != n-2) {PetscObjectDereference((PetscObject)dB);}
631: dB = B;
632: }
633: if (n > 1) {PetscObjectDereference((PetscObject)dB);}
634: } else {
635: for (i=n-2; i>-1; i--) {
636: KSPGetOperators(mglevels[i]->smoothd,NULL,&B,NULL);
637: MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);
638: KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);
639: dB = B;
640: }
641: }
642: } else if (!mg->galerkin && pc->dm && pc->dm->x) {
643: /* need to restrict Jacobian location to coarser meshes for evaluation */
644: for (i=n-2; i>-1; i--) {
645: Mat R;
646: Vec rscale;
647: if (!mglevels[i]->smoothd->dm->x) {
648: Vec *vecs;
649: KSPGetVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);
651: mglevels[i]->smoothd->dm->x = vecs[0];
653: PetscFree(vecs);
654: }
655: PCMGGetRestriction(pc,i+1,&R);
656: PCMGGetRScale(pc,i+1,&rscale);
657: MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);
658: VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);
659: }
660: }
661: if (!mg->galerkin && pc->dm) {
662: for (i=n-2; i>=0; i--) {
663: DM dmfine,dmcoarse;
664: Mat Restrict,Inject;
665: Vec rscale;
666: KSPGetDM(mglevels[i+1]->smoothd,&dmfine);
667: KSPGetDM(mglevels[i]->smoothd,&dmcoarse);
668: PCMGGetRestriction(pc,i+1,&Restrict);
669: PCMGGetRScale(pc,i+1,&rscale);
670: Inject = NULL; /* Callback should create it if it needs Injection */
671: DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);
672: }
673: }
675: if (!pc->setupcalled) {
676: for (i=0; i<n; i++) {
677: KSPSetFromOptions(mglevels[i]->smoothd);
678: }
679: for (i=1; i<n; i++) {
680: if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
681: KSPSetFromOptions(mglevels[i]->smoothu);
682: }
683: }
684: for (i=1; i<n; i++) {
685: PCMGGetInterpolation(pc,i,&mglevels[i]->interpolate);
686: PCMGGetRestriction(pc,i,&mglevels[i]->restrct);
687: }
688: for (i=0; i<n-1; i++) {
689: if (!mglevels[i]->b) {
690: Vec *vec;
691: KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,NULL);
692: PCMGSetRhs(pc,i,*vec);
693: VecDestroy(vec);
694: PetscFree(vec);
695: }
696: if (!mglevels[i]->r && i) {
697: VecDuplicate(mglevels[i]->b,&tvec);
698: PCMGSetR(pc,i,tvec);
699: VecDestroy(&tvec);
700: }
701: if (!mglevels[i]->x) {
702: VecDuplicate(mglevels[i]->b,&tvec);
703: PCMGSetX(pc,i,tvec);
704: VecDestroy(&tvec);
705: }
706: }
707: if (n != 1 && !mglevels[n-1]->r) {
708: /* PCMGSetR() on the finest level if user did not supply it */
709: Vec *vec;
710: KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);
711: PCMGSetR(pc,n-1,*vec);
712: VecDestroy(vec);
713: PetscFree(vec);
714: }
715: }
717: if (pc->dm) {
718: /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
719: for (i=0; i<n-1; i++) {
720: if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
721: }
722: }
724: for (i=1; i<n; i++) {
725: if (mglevels[i]->smoothu == mglevels[i]->smoothd) {
726: /* if doing only down then initial guess is zero */
727: KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
728: }
729: if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
730: KSPSetUp(mglevels[i]->smoothd);
731: if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
732: if (!mglevels[i]->residual) {
733: Mat mat;
734: KSPGetOperators(mglevels[i]->smoothd,NULL,&mat,NULL);
735: PCMGSetResidual(pc,i,PCMGResidual_Default,mat);
736: }
737: }
738: for (i=1; i<n; i++) {
739: if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
740: Mat downmat,downpmat;
741: MatStructure matflag;
743: /* check if operators have been set for up, if not use down operators to set them */
744: KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);
745: if (!opsset) {
746: KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);
747: KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);
748: }
750: KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);
751: if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
752: KSPSetUp(mglevels[i]->smoothu);
753: if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
754: }
755: }
757: /*
758: If coarse solver is not direct method then DO NOT USE preonly
759: */
760: PetscObjectTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);
761: if (preonly) {
762: PetscObjectTypeCompare((PetscObject)cpc,PCLU,&lu);
763: PetscObjectTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);
764: PetscObjectTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);
765: PetscObjectTypeCompare((PetscObject)cpc,PCSVD,&svd);
766: if (!lu && !redundant && !cholesky && !svd) {
767: KSPSetType(mglevels[0]->smoothd,KSPGMRES);
768: }
769: }
771: if (!pc->setupcalled) {
772: KSPSetFromOptions(mglevels[0]->smoothd);
773: }
775: if (mglevels[0]->eventsmoothsetup) {PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);}
776: KSPSetUp(mglevels[0]->smoothd);
777: if (mglevels[0]->eventsmoothsetup) {PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);}
779: /*
780: Dump the interpolation/restriction matrices plus the
781: Jacobian/stiffness on each level. This allows MATLAB users to
782: easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
784: Only support one or the other at the same time.
785: */
786: #if defined(PETSC_USE_SOCKET_VIEWER)
787: PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);
788: if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
789: dump = PETSC_FALSE;
790: #endif
791: PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);
792: if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
794: if (viewer) {
795: for (i=1; i<n; i++) {
796: MatView(mglevels[i]->restrct,viewer);
797: }
798: for (i=0; i<n; i++) {
799: KSPGetPC(mglevels[i]->smoothd,&pc);
800: MatView(pc->mat,viewer);
801: }
802: }
803: return(0);
804: }
806: /* -------------------------------------------------------------------------------------*/
810: /*@
811: PCMGGetLevels - Gets the number of levels to use with MG.
813: Not Collective
815: Input Parameter:
816: . pc - the preconditioner context
818: Output parameter:
819: . levels - the number of levels
821: Level: advanced
823: .keywords: MG, get, levels, multigrid
825: .seealso: PCMGSetLevels()
826: @*/
827: PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
828: {
829: PC_MG *mg = (PC_MG*)pc->data;
834: *levels = mg->nlevels;
835: return(0);
836: }
840: /*@
841: PCMGSetType - Determines the form of multigrid to use:
842: multiplicative, additive, full, or the Kaskade algorithm.
844: Logically Collective on PC
846: Input Parameters:
847: + pc - the preconditioner context
848: - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
849: PC_MG_FULL, PC_MG_KASKADE
851: Options Database Key:
852: . -pc_mg_type <form> - Sets <form>, one of multiplicative,
853: additive, full, kaskade
855: Level: advanced
857: .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid
859: .seealso: PCMGSetLevels()
860: @*/
861: PetscErrorCode PCMGSetType(PC pc,PCMGType form)
862: {
863: PC_MG *mg = (PC_MG*)pc->data;
868: mg->am = form;
869: if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
870: else pc->ops->applyrichardson = 0;
871: return(0);
872: }
876: /*@
877: PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more
878: complicated cycling.
880: Logically Collective on PC
882: Input Parameters:
883: + pc - the multigrid context
884: - PC_MG_CYCLE_V or PC_MG_CYCLE_W
886: Options Database Key:
887: $ -pc_mg_cycle_type v or w
889: Level: advanced
891: .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
893: .seealso: PCMGSetCycleTypeOnLevel()
894: @*/
895: PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n)
896: {
897: PC_MG *mg = (PC_MG*)pc->data;
898: PC_MG_Levels **mglevels = mg->levels;
899: PetscInt i,levels;
903: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
905: levels = mglevels[0]->levels;
907: for (i=0; i<levels; i++) mglevels[i]->cycles = n;
908: return(0);
909: }
913: /*@
914: PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
915: of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
917: Logically Collective on PC
919: Input Parameters:
920: + pc - the multigrid context
921: - n - number of cycles (default is 1)
923: Options Database Key:
924: $ -pc_mg_multiplicative_cycles n
926: Level: advanced
928: Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
930: .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
932: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
933: @*/
934: PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
935: {
936: PC_MG *mg = (PC_MG*)pc->data;
937: PC_MG_Levels **mglevels = mg->levels;
938: PetscInt i,levels;
942: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
944: levels = mglevels[0]->levels;
946: for (i=0; i<levels; i++) mg->cyclesperpcapply = n;
947: return(0);
948: }
952: /*@
953: PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
954: finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t
956: Logically Collective on PC
958: Input Parameters:
959: + pc - the multigrid context
960: - use - PETSC_TRUE to use the Galerkin process to compute coarse-level operators
962: Options Database Key:
963: $ -pc_mg_galerkin
965: Level: intermediate
967: .keywords: MG, set, Galerkin
969: .seealso: PCMGGetGalerkin()
971: @*/
972: PetscErrorCode PCMGSetGalerkin(PC pc,PetscBool use)
973: {
974: PC_MG *mg = (PC_MG*)pc->data;
978: mg->galerkin = use ? 1 : 0;
979: return(0);
980: }
984: /*@
985: PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
986: A_i-1 = r_i * A_i * r_i^t
988: Not Collective
990: Input Parameter:
991: . pc - the multigrid context
993: Output Parameter:
994: . gelerkin - PETSC_TRUE or PETSC_FALSE
996: Options Database Key:
997: $ -pc_mg_galerkin
999: Level: intermediate
1001: .keywords: MG, set, Galerkin
1003: .seealso: PCMGSetGalerkin()
1005: @*/
1006: PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin)
1007: {
1008: PC_MG *mg = (PC_MG*)pc->data;
1012: *galerkin = (PetscBool)mg->galerkin;
1013: return(0);
1014: }
1018: /*@
1019: PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to
1020: use on all levels. Use PCMGGetSmootherDown() to set different
1021: pre-smoothing steps on different levels.
1023: Logically Collective on PC
1025: Input Parameters:
1026: + mg - the multigrid context
1027: - n - the number of smoothing steps
1029: Options Database Key:
1030: . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps
1032: Level: advanced
1034: .keywords: MG, smooth, down, pre-smoothing, steps, multigrid
1036: .seealso: PCMGSetNumberSmoothUp()
1037: @*/
1038: PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n)
1039: {
1040: PC_MG *mg = (PC_MG*)pc->data;
1041: PC_MG_Levels **mglevels = mg->levels;
1043: PetscInt i,levels;
1047: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1049: levels = mglevels[0]->levels;
1051: for (i=1; i<levels; i++) {
1052: /* make sure smoother up and down are different */
1053: PCMGGetSmootherUp(pc,i,NULL);
1054: KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1056: mg->default_smoothd = n;
1057: }
1058: return(0);
1059: }
1063: /*@
1064: PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use
1065: on all levels. Use PCMGGetSmootherUp() to set different numbers of
1066: post-smoothing steps on different levels.
1068: Logically Collective on PC
1070: Input Parameters:
1071: + mg - the multigrid context
1072: - n - the number of smoothing steps
1074: Options Database Key:
1075: . -pc_mg_smoothup <n> - Sets number of post-smoothing steps
1077: Level: advanced
1079: Note: this does not set a value on the coarsest grid, since we assume that
1080: there is no separate smooth up on the coarsest grid.
1082: .keywords: MG, smooth, up, post-smoothing, steps, multigrid
1084: .seealso: PCMGSetNumberSmoothDown()
1085: @*/
1086: PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n)
1087: {
1088: PC_MG *mg = (PC_MG*)pc->data;
1089: PC_MG_Levels **mglevels = mg->levels;
1091: PetscInt i,levels;
1095: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1097: levels = mglevels[0]->levels;
1099: for (i=1; i<levels; i++) {
1100: /* make sure smoother up and down are different */
1101: PCMGGetSmootherUp(pc,i,NULL);
1102: KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1104: mg->default_smoothu = n;
1105: }
1106: return(0);
1107: }
1109: /* ----------------------------------------------------------------------------------------*/
1111: /*MC
1112: PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1113: information about the coarser grid matrices and restriction/interpolation operators.
1115: Options Database Keys:
1116: + -pc_mg_levels <nlevels> - number of levels including finest
1117: . -pc_mg_cycles v or w
1118: . -pc_mg_smoothup <n> - number of smoothing steps after interpolation
1119: . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator
1120: . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1121: . -pc_mg_log - log information about time spent on each level of the solver
1122: . -pc_mg_monitor - print information on the multigrid convergence
1123: . -pc_mg_galerkin - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1124: . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1125: . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1126: to the Socket viewer for reading from MATLAB.
1127: - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1128: to the binary output file called binaryoutput
1130: Notes: By default this uses GMRES on the fine grid smoother so this should be used with KSPFGMRES or the smoother changed to not use GMRES
1132: When run with a single level the smoother options are used on that level NOT the coarse grid solver options
1134: Level: intermediate
1136: Concepts: multigrid/multilevel
1138: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1139: PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(),
1140: PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1141: PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1142: PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1143: M*/
1147: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1148: {
1149: PC_MG *mg;
1153: PetscNewLog(pc,PC_MG,&mg);
1154: pc->data = (void*)mg;
1155: mg->nlevels = -1;
1157: pc->useAmat = PETSC_TRUE;
1159: pc->ops->apply = PCApply_MG;
1160: pc->ops->setup = PCSetUp_MG;
1161: pc->ops->reset = PCReset_MG;
1162: pc->ops->destroy = PCDestroy_MG;
1163: pc->ops->setfromoptions = PCSetFromOptions_MG;
1164: pc->ops->view = PCView_MG;
1165: return(0);
1166: }