Actual source code: mg.c
petsc-3.15.0 2021-03-30
2: /*
3: Defines the multigrid preconditioner interface.
4: */
5: #include <petsc/private/pcmgimpl.h>
6: #include <petsc/private/kspimpl.h>
7: #include <petscdm.h>
8: PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*);
10: /*
11: Contains the list of registered coarse space construction routines
12: */
13: PetscFunctionList PCMGCoarseList = NULL;
15: PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PetscBool transpose,PetscBool matapp,PCRichardsonConvergedReason *reason)
16: {
17: PC_MG *mg = (PC_MG*)pc->data;
18: PC_MG_Levels *mgc,*mglevels = *mglevelsin;
20: PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles;
23: if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
24: if (!transpose) {
25: if (matapp) {
26: KSPMatSolve(mglevels->smoothd,mglevels->B,mglevels->X); /* pre-smooth */
27: KSPCheckSolve(mglevels->smoothd,pc,NULL);
28: } else {
29: KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x); /* pre-smooth */
30: KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);
31: }
32: } else {
33: if (matapp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported");
34: KSPSolveTranspose(mglevels->smoothu,mglevels->b,mglevels->x); /* transpose of post-smooth */
35: KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);
36: }
37: if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
38: if (mglevels->level) { /* not the coarsest grid */
39: if (mglevels->eventresidual) {PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);}
40: if (matapp && !mglevels->R) {
41: MatDuplicate(mglevels->B,MAT_DO_NOT_COPY_VALUES,&mglevels->R);
42: }
43: if (!transpose) {
44: if (matapp) { (*mglevels->matresidual)(mglevels->A,mglevels->B,mglevels->X,mglevels->R); }
45: else { (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r); }
46: } else {
47: if (matapp) { (*mglevels->matresidualtranspose)(mglevels->A,mglevels->B,mglevels->X,mglevels->R); }
48: else { (*mglevels->residualtranspose)(mglevels->A,mglevels->b,mglevels->x,mglevels->r); }
49: }
50: if (mglevels->eventresidual) {PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);}
52: /* if on finest level and have convergence criteria set */
53: if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) {
54: PetscReal rnorm;
55: VecNorm(mglevels->r,NORM_2,&rnorm);
56: if (rnorm <= mg->ttol) {
57: if (rnorm < mg->abstol) {
58: *reason = PCRICHARDSON_CONVERGED_ATOL;
59: PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);
60: } else {
61: *reason = PCRICHARDSON_CONVERGED_RTOL;
62: 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);
63: }
64: return(0);
65: }
66: }
68: mgc = *(mglevelsin - 1);
69: if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
70: if (!transpose) {
71: if (matapp) { MatMatRestrict(mglevels->restrct,mglevels->R,&mgc->B); }
72: else { MatRestrict(mglevels->restrct,mglevels->r,mgc->b); }
73: } else {
74: if (matapp) { MatMatRestrict(mglevels->interpolate,mglevels->R,&mgc->B); }
75: else { MatRestrict(mglevels->interpolate,mglevels->r,mgc->b); }
76: }
77: if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
78: if (matapp) {
79: if (!mgc->X) {
80: MatDuplicate(mgc->B,MAT_DO_NOT_COPY_VALUES,&mgc->X);
81: } else {
82: MatZeroEntries(mgc->X);
83: }
84: } else {
85: VecZeroEntries(mgc->x);
86: }
87: while (cycles--) {
88: PCMGMCycle_Private(pc,mglevelsin-1,transpose,matapp,reason);
89: }
90: if (mglevels->eventinterprestrict) {PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);}
91: if (!transpose) {
92: if (matapp) { MatMatInterpolateAdd(mglevels->interpolate,mgc->X,mglevels->X,&mglevels->X); }
93: else { MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x); }
94: } else {
95: MatInterpolateAdd(mglevels->restrct,mgc->x,mglevels->x,mglevels->x);
96: }
97: if (mglevels->eventinterprestrict) {PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);}
98: if (mglevels->eventsmoothsolve) {PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);}
99: if (!transpose) {
100: if (matapp) {
101: KSPMatSolve(mglevels->smoothu,mglevels->B,mglevels->X); /* post smooth */
102: KSPCheckSolve(mglevels->smoothu,pc,NULL);
103: } else {
104: KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x); /* post smooth */
105: KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);
106: }
107: } else {
108: if (matapp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported");
109: KSPSolveTranspose(mglevels->smoothd,mglevels->b,mglevels->x); /* post smooth */
110: KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);
111: }
112: if (mglevels->cr) {
113: if (matapp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported");
114: /* TODO Turn on copy and turn off noisy if we have an exact solution
115: VecCopy(mglevels->x, mglevels->crx);
116: VecCopy(mglevels->b, mglevels->crb); */
117: KSPSetNoisy_Private(mglevels->crx);
118: KSPSolve(mglevels->cr,mglevels->crb,mglevels->crx); /* compatible relaxation */
119: KSPCheckSolve(mglevels->cr,pc,mglevels->crx);
120: }
121: if (mglevels->eventsmoothsolve) {PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);}
122: }
123: return(0);
124: }
126: 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)
127: {
128: PC_MG *mg = (PC_MG*)pc->data;
129: PC_MG_Levels **mglevels = mg->levels;
131: PC tpc;
132: PetscBool changeu,changed;
133: PetscInt levels = mglevels[0]->levels,i;
136: /* When the DM is supplying the matrix then it will not exist until here */
137: for (i=0; i<levels; i++) {
138: if (!mglevels[i]->A) {
139: KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
140: PetscObjectReference((PetscObject)mglevels[i]->A);
141: }
142: }
144: KSPGetPC(mglevels[levels-1]->smoothd,&tpc);
145: PCPreSolveChangeRHS(tpc,&changed);
146: KSPGetPC(mglevels[levels-1]->smoothu,&tpc);
147: PCPreSolveChangeRHS(tpc,&changeu);
148: if (!changed && !changeu) {
149: VecDestroy(&mglevels[levels-1]->b);
150: mglevels[levels-1]->b = b;
151: } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
152: if (!mglevels[levels-1]->b) {
153: Vec *vec;
155: KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);
156: mglevels[levels-1]->b = *vec;
157: PetscFree(vec);
158: }
159: VecCopy(b,mglevels[levels-1]->b);
160: }
161: mglevels[levels-1]->x = x;
163: mg->rtol = rtol;
164: mg->abstol = abstol;
165: mg->dtol = dtol;
166: if (rtol) {
167: /* compute initial residual norm for relative convergence test */
168: PetscReal rnorm;
169: if (zeroguess) {
170: VecNorm(b,NORM_2,&rnorm);
171: } else {
172: (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);
173: VecNorm(w,NORM_2,&rnorm);
174: }
175: mg->ttol = PetscMax(rtol*rnorm,abstol);
176: } else if (abstol) mg->ttol = abstol;
177: else mg->ttol = 0.0;
179: /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
180: stop prematurely due to small residual */
181: for (i=1; i<levels; i++) {
182: KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
183: if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
184: /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
185: KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
186: KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
187: }
188: }
190: *reason = (PCRichardsonConvergedReason)0;
191: for (i=0; i<its; i++) {
192: PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_FALSE,PETSC_FALSE,reason);
193: if (*reason) break;
194: }
195: if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS;
196: *outits = i;
197: if (!changed && !changeu) mglevels[levels-1]->b = NULL;
198: return(0);
199: }
201: PetscErrorCode PCReset_MG(PC pc)
202: {
203: PC_MG *mg = (PC_MG*)pc->data;
204: PC_MG_Levels **mglevels = mg->levels;
206: PetscInt i,c,n;
209: if (mglevels) {
210: n = mglevels[0]->levels;
211: for (i=0; i<n-1; i++) {
212: VecDestroy(&mglevels[i+1]->r);
213: VecDestroy(&mglevels[i]->b);
214: VecDestroy(&mglevels[i]->x);
215: MatDestroy(&mglevels[i+1]->R);
216: MatDestroy(&mglevels[i]->B);
217: MatDestroy(&mglevels[i]->X);
218: VecDestroy(&mglevels[i]->crx);
219: VecDestroy(&mglevels[i]->crb);
220: MatDestroy(&mglevels[i+1]->restrct);
221: MatDestroy(&mglevels[i+1]->interpolate);
222: MatDestroy(&mglevels[i+1]->inject);
223: VecDestroy(&mglevels[i+1]->rscale);
224: }
225: VecDestroy(&mglevels[n-1]->crx);
226: VecDestroy(&mglevels[n-1]->crb);
227: /* this is not null only if the smoother on the finest level
228: changes the rhs during PreSolve */
229: VecDestroy(&mglevels[n-1]->b);
230: MatDestroy(&mglevels[n-1]->B);
232: for (i=0; i<n; i++) {
233: if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {VecDestroy(&mglevels[i]->coarseSpace[c]);}
234: PetscFree(mglevels[i]->coarseSpace);
235: mglevels[i]->coarseSpace = NULL;
236: MatDestroy(&mglevels[i]->A);
237: if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
238: KSPReset(mglevels[i]->smoothd);
239: }
240: KSPReset(mglevels[i]->smoothu);
241: if (mglevels[i]->cr) {KSPReset(mglevels[i]->cr);}
242: }
243: mg->Nc = 0;
244: }
245: return(0);
246: }
248: /* Implementing CR
250: We only want to make corrections that ``do not change'' the coarse solution. What we mean by not changing is that if I prolong my coarse solution to the fine grid and then inject that fine solution back to the coarse grid, I get the same answer. Injection is what Brannick calls R. We want the complementary projector to Inj, which we will call S, after Brannick, so that Inj S = 0. Now the orthogonal projector onto the range of Inj^T is
252: Inj^T (Inj Inj^T)^{-1} Inj
254: and if Inj is a VecScatter, as it is now in PETSc, we have
256: Inj^T Inj
258: and
260: S = I - Inj^T Inj
262: since
264: Inj S = Inj - (Inj Inj^T) Inj = 0.
266: Brannick suggests
268: A \to S^T A S \qquad\mathrm{and}\qquad M \to S^T M S
270: but I do not think his :math:`S^T S = I` is correct. Our S is an orthogonal projector, so :math:`S^T S = S^2 = S`. We will use
272: M^{-1} A \to S M^{-1} A S
274: In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left.
276: Check: || Inj P - I ||_F < tol
277: Check: In general, Inj Inj^T = I
278: */
280: typedef struct {
281: PC mg; /* The PCMG object */
282: PetscInt l; /* The multigrid level for this solver */
283: Mat Inj; /* The injection matrix */
284: Mat S; /* I - Inj^T Inj */
285: } CRContext;
287: static PetscErrorCode CRSetup_Private(PC pc)
288: {
289: CRContext *ctx;
290: Mat It;
294: PCShellGetContext(pc, (void **) &ctx);
295: PCMGGetInjection(ctx->mg, ctx->l, &It);
296: if (!It) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG");
297: MatCreateTranspose(It, &ctx->Inj);
298: MatCreateNormal(ctx->Inj, &ctx->S);
299: MatScale(ctx->S, -1.0);
300: MatShift(ctx->S, 1.0);
301: return(0);
302: }
304: static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y)
305: {
306: CRContext *ctx;
310: PCShellGetContext(pc, (void **) &ctx);
311: MatMult(ctx->S, x, y);
312: return(0);
313: }
315: static PetscErrorCode CRDestroy_Private(PC pc)
316: {
317: CRContext *ctx;
321: PCShellGetContext(pc, (void **) &ctx);
322: MatDestroy(&ctx->Inj);
323: MatDestroy(&ctx->S);
324: PetscFree(ctx);
325: PCShellSetContext(pc, NULL);
326: return(0);
327: }
329: static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr)
330: {
331: CRContext *ctx;
335: PCCreate(PetscObjectComm((PetscObject) pc), cr);
336: PetscObjectSetName((PetscObject) *cr, "S (complementary projector to injection)");
337: PetscCalloc1(1, &ctx);
338: ctx->mg = pc;
339: ctx->l = l;
340: PCSetType(*cr, PCSHELL);
341: PCShellSetContext(*cr, ctx);
342: PCShellSetApply(*cr, CRApply_Private);
343: PCShellSetSetUp(*cr, CRSetup_Private);
344: PCShellSetDestroy(*cr, CRDestroy_Private);
345: return(0);
346: }
348: PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms)
349: {
351: PC_MG *mg = (PC_MG*)pc->data;
352: MPI_Comm comm;
353: PC_MG_Levels **mglevels = mg->levels;
354: PCMGType mgtype = mg->am;
355: PetscInt mgctype = (PetscInt) PC_MG_CYCLE_V;
356: PetscInt i;
357: PetscMPIInt size;
358: const char *prefix;
359: PC ipc;
360: PetscInt n;
365: if (mg->nlevels == levels) return(0);
366: PetscObjectGetComm((PetscObject)pc,&comm);
367: if (mglevels) {
368: mgctype = mglevels[0]->cycles;
369: /* changing the number of levels so free up the previous stuff */
370: PCReset_MG(pc);
371: n = mglevels[0]->levels;
372: for (i=0; i<n; i++) {
373: if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
374: KSPDestroy(&mglevels[i]->smoothd);
375: }
376: KSPDestroy(&mglevels[i]->smoothu);
377: KSPDestroy(&mglevels[i]->cr);
378: PetscFree(mglevels[i]);
379: }
380: PetscFree(mg->levels);
381: }
383: mg->nlevels = levels;
385: PetscMalloc1(levels,&mglevels);
386: PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));
388: PCGetOptionsPrefix(pc,&prefix);
390: mg->stageApply = 0;
391: for (i=0; i<levels; i++) {
392: PetscNewLog(pc,&mglevels[i]);
394: mglevels[i]->level = i;
395: mglevels[i]->levels = levels;
396: mglevels[i]->cycles = mgctype;
397: mg->default_smoothu = 2;
398: mg->default_smoothd = 2;
399: mglevels[i]->eventsmoothsetup = 0;
400: mglevels[i]->eventsmoothsolve = 0;
401: mglevels[i]->eventresidual = 0;
402: mglevels[i]->eventinterprestrict = 0;
404: if (comms) comm = comms[i];
405: if (comm != MPI_COMM_NULL) {
406: KSPCreate(comm,&mglevels[i]->smoothd);
407: KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);
408: PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);
409: KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);
410: PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);
411: if (i || levels == 1) {
412: char tprefix[128];
414: KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);
415: KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);
416: KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);
417: KSPGetPC(mglevels[i]->smoothd,&ipc);
418: PCSetType(ipc,PCSOR);
419: KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);
421: PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i);
422: KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);
423: } else {
424: KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");
426: /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
427: KSPSetType(mglevels[0]->smoothd,KSPPREONLY);
428: KSPGetPC(mglevels[0]->smoothd,&ipc);
429: MPI_Comm_size(comm,&size);
430: if (size > 1) {
431: PCSetType(ipc,PCREDUNDANT);
432: } else {
433: PCSetType(ipc,PCLU);
434: }
435: PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);
436: }
437: PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);
438: }
439: mglevels[i]->smoothu = mglevels[i]->smoothd;
440: mg->rtol = 0.0;
441: mg->abstol = 0.0;
442: mg->dtol = 0.0;
443: mg->ttol = 0.0;
444: mg->cyclesperpcapply = 1;
445: }
446: mg->levels = mglevels;
447: PCMGSetType(pc,mgtype);
448: return(0);
449: }
451: /*@C
452: PCMGSetLevels - Sets the number of levels to use with MG.
453: Must be called before any other MG routine.
455: Logically Collective on PC
457: Input Parameters:
458: + pc - the preconditioner context
459: . levels - the number of levels
460: - comms - optional communicators for each level; this is to allow solving the coarser problems
461: on smaller sets of processes. For processes that are not included in the computation
462: you must pass MPI_COMM_NULL.
464: Level: intermediate
466: Notes:
467: If the number of levels is one then the multigrid uses the -mg_levels prefix
468: for setting the level options rather than the -mg_coarse prefix.
470: You can free the information in comms after this routine is called.
472: The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level
473: are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on
474: the two levels, rank 0 in the original communicator will pass in an array of 2 communicators
475: of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators
476: the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate
477: in the coarse grid solve.
479: Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one
480: must take special care in providing the restriction and interpolation operation. We recommend
481: providing these as two step operations; first perform a standard restriction or interpolation on
482: the full number of ranks for that level and then use an MPI call to copy the resulting vector
483: array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, not in both
484: cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and
485: recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries.
488: .seealso: PCMGSetType(), PCMGGetLevels()
489: @*/
490: PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms)
491: {
497: PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));
498: return(0);
499: }
502: PetscErrorCode PCDestroy_MG(PC pc)
503: {
505: PC_MG *mg = (PC_MG*)pc->data;
506: PC_MG_Levels **mglevels = mg->levels;
507: PetscInt i,n;
510: PCReset_MG(pc);
511: if (mglevels) {
512: n = mglevels[0]->levels;
513: for (i=0; i<n; i++) {
514: if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
515: KSPDestroy(&mglevels[i]->smoothd);
516: }
517: KSPDestroy(&mglevels[i]->smoothu);
518: KSPDestroy(&mglevels[i]->cr);
519: PetscFree(mglevels[i]);
520: }
521: PetscFree(mg->levels);
522: }
523: PetscFree(pc->data);
524: PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);
525: PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);
526: return(0);
527: }
530: /*
531: PCApply_MG - Runs either an additive, multiplicative, Kaskadic
532: or full cycle of multigrid.
534: Note:
535: A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
536: */
537: static PetscErrorCode PCApply_MG_Internal(PC pc,Vec b,Vec x,Mat B,Mat X,PetscBool transpose)
538: {
539: PC_MG *mg = (PC_MG*)pc->data;
540: PC_MG_Levels **mglevels = mg->levels;
542: PC tpc;
543: PetscInt levels = mglevels[0]->levels,i;
544: PetscBool changeu,changed,matapp;
547: matapp = (PetscBool)(B && X);
548: if (mg->stageApply) {PetscLogStagePush(mg->stageApply);}
549: /* When the DM is supplying the matrix then it will not exist until here */
550: for (i=0; i<levels; i++) {
551: if (!mglevels[i]->A) {
552: KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);
553: PetscObjectReference((PetscObject)mglevels[i]->A);
554: }
555: }
557: KSPGetPC(mglevels[levels-1]->smoothd,&tpc);
558: PCPreSolveChangeRHS(tpc,&changed);
559: KSPGetPC(mglevels[levels-1]->smoothu,&tpc);
560: PCPreSolveChangeRHS(tpc,&changeu);
561: if (!changeu && !changed) {
562: if (matapp) {
563: MatDestroy(&mglevels[levels-1]->B);
564: mglevels[levels-1]->B = B;
565: } else {
566: VecDestroy(&mglevels[levels-1]->b);
567: mglevels[levels-1]->b = b;
568: }
569: } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
570: if (matapp) {
571: if (mglevels[levels-1]->B) {
572: PetscInt N1,N2;
573: PetscBool flg;
575: MatGetSize(mglevels[levels-1]->B,NULL,&N1);
576: MatGetSize(B,NULL,&N2);
577: PetscObjectTypeCompare((PetscObject)mglevels[levels-1]->B,((PetscObject)B)->type_name,&flg);
578: if (N1 != N2 || !flg) {
579: MatDestroy(&mglevels[levels-1]->B);
580: }
581: }
582: if (!mglevels[levels-1]->B) {
583: MatDuplicate(B,MAT_COPY_VALUES,&mglevels[levels-1]->B);
584: } else {
585: MatCopy(B,mglevels[levels-1]->B,SAME_NONZERO_PATTERN);
586: }
587: } else {
588: if (!mglevels[levels-1]->b) {
589: Vec *vec;
591: KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);
592: mglevels[levels-1]->b = *vec;
593: PetscFree(vec);
594: }
595: VecCopy(b,mglevels[levels-1]->b);
596: }
597: }
598: if (matapp) { mglevels[levels-1]->X = X; }
599: else { mglevels[levels-1]->x = x; }
601: /* If coarser Xs are present, it means we have already block applied the PC at least once
602: Reset operators if sizes/type do no match */
603: if (matapp && levels > 1 && mglevels[levels-2]->X) {
604: PetscInt Xc,Bc;
605: PetscBool flg;
607: MatGetSize(mglevels[levels-2]->X,NULL,&Xc);
608: MatGetSize(mglevels[levels-1]->B,NULL,&Bc);
609: PetscObjectTypeCompare((PetscObject)mglevels[levels-2]->X,((PetscObject)mglevels[levels-1]->X)->type_name,&flg);
610: if (Xc != Bc || !flg) {
611: MatDestroy(&mglevels[levels-1]->R);
612: for (i=0;i<levels-1;i++) {
613: MatDestroy(&mglevels[i]->R);
614: MatDestroy(&mglevels[i]->B);
615: MatDestroy(&mglevels[i]->X);
616: }
617: }
618: }
620: if (mg->am == PC_MG_MULTIPLICATIVE) {
621: if (matapp) { MatZeroEntries(X); }
622: else { VecZeroEntries(x); }
623: for (i=0; i<mg->cyclesperpcapply; i++) {
624: PCMGMCycle_Private(pc,mglevels+levels-1,transpose,matapp,NULL);
625: }
626: } else if (mg->am == PC_MG_ADDITIVE) {
627: PCMGACycle_Private(pc,mglevels,transpose,matapp);
628: } else if (mg->am == PC_MG_KASKADE) {
629: PCMGKCycle_Private(pc,mglevels,transpose,matapp);
630: } else {
631: PCMGFCycle_Private(pc,mglevels,transpose,matapp);
632: }
633: if (mg->stageApply) {PetscLogStagePop();}
634: if (!changeu && !changed) {
635: if (matapp) { mglevels[levels-1]->B = NULL; }
636: else { mglevels[levels-1]->b = NULL; }
637: }
638: return(0);
639: }
641: static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x)
642: {
646: PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_FALSE);
647: return(0);
648: }
650: static PetscErrorCode PCApplyTranspose_MG(PC pc,Vec b,Vec x)
651: {
655: PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_TRUE);
656: return(0);
657: }
659: static PetscErrorCode PCMatApply_MG(PC pc,Mat b,Mat x)
660: {
664: PCApply_MG_Internal(pc,NULL,NULL,b,x,PETSC_FALSE);
665: return(0);
666: }
668: PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc)
669: {
670: PetscErrorCode ierr;
671: PetscInt levels,cycles;
672: PetscBool flg, flg2;
673: PC_MG *mg = (PC_MG*)pc->data;
674: PC_MG_Levels **mglevels;
675: PCMGType mgtype;
676: PCMGCycleType mgctype;
677: PCMGGalerkinType gtype;
680: levels = PetscMax(mg->nlevels,1);
681: PetscOptionsHead(PetscOptionsObject,"Multigrid options");
682: PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);
683: if (!flg && !mg->levels && pc->dm) {
684: DMGetRefineLevel(pc->dm,&levels);
685: levels++;
686: mg->usedmfornumberoflevels = PETSC_TRUE;
687: }
688: PCMGSetLevels(pc,levels,NULL);
689: mglevels = mg->levels;
691: mgctype = (PCMGCycleType) mglevels[0]->cycles;
692: PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);
693: if (flg) {
694: PCMGSetCycleType(pc,mgctype);
695: }
696: gtype = mg->galerkin;
697: PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)>ype,&flg);
698: if (flg) {
699: PCMGSetGalerkin(pc,gtype);
700: }
701: flg2 = PETSC_FALSE;
702: PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);
703: if (flg) {PCMGSetAdaptInterpolation(pc, flg2);}
704: PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg);
705: 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);
706: PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);
707: flg2 = PETSC_FALSE;
708: PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg);
709: if (flg) {PCMGSetAdaptCR(pc, flg2);}
710: flg = PETSC_FALSE;
711: PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);
712: if (flg) {
713: PCMGSetDistinctSmoothUp(pc);
714: }
715: mgtype = mg->am;
716: PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);
717: if (flg) {
718: PCMGSetType(pc,mgtype);
719: }
720: if (mg->am == PC_MG_MULTIPLICATIVE) {
721: PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);
722: if (flg) {
723: PCMGMultiplicativeSetCycles(pc,cycles);
724: }
725: }
726: flg = PETSC_FALSE;
727: PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);
728: if (flg) {
729: PetscInt i;
730: char eventname[128];
732: levels = mglevels[0]->levels;
733: for (i=0; i<levels; i++) {
734: sprintf(eventname,"MGSetup Level %d",(int)i);
735: PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);
736: sprintf(eventname,"MGSmooth Level %d",(int)i);
737: PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);
738: if (i) {
739: sprintf(eventname,"MGResid Level %d",(int)i);
740: PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);
741: sprintf(eventname,"MGInterp Level %d",(int)i);
742: PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);
743: }
744: }
746: #if defined(PETSC_USE_LOG)
747: {
748: const char *sname = "MG Apply";
749: PetscStageLog stageLog;
750: PetscInt st;
752: PetscLogGetStageLog(&stageLog);
753: for (st = 0; st < stageLog->numStages; ++st) {
754: PetscBool same;
756: PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);
757: if (same) mg->stageApply = st;
758: }
759: if (!mg->stageApply) {
760: PetscLogStageRegister(sname, &mg->stageApply);
761: }
762: }
763: #endif
764: }
765: PetscOptionsTail();
766: /* Check option consistency */
767: PCMGGetGalerkin(pc, >ype);
768: PCMGGetAdaptInterpolation(pc, &flg);
769: if (flg && (gtype >= PC_MG_GALERKIN_NONE)) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
770: return(0);
771: }
773: const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL};
774: const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL};
775: const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL};
776: const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL};
778: #include <petscdraw.h>
779: PetscErrorCode PCView_MG(PC pc,PetscViewer viewer)
780: {
781: PC_MG *mg = (PC_MG*)pc->data;
782: PC_MG_Levels **mglevels = mg->levels;
784: PetscInt levels = mglevels ? mglevels[0]->levels : 0,i;
785: PetscBool iascii,isbinary,isdraw;
788: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
789: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);
790: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);
791: if (iascii) {
792: const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";
793: PetscViewerASCIIPrintf(viewer," type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);
794: if (mg->am == PC_MG_MULTIPLICATIVE) {
795: PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);
796: }
797: if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
798: PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");
799: } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
800: PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for pmat\n");
801: } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
802: PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for mat\n");
803: } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
804: PetscViewerASCIIPrintf(viewer," Using externally compute Galerkin coarse grid matrices\n");
805: } else {
806: PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");
807: }
808: if (mg->view){
809: (*mg->view)(pc,viewer);
810: }
811: for (i=0; i<levels; i++) {
812: if (!i) {
813: PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);
814: } else {
815: PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);
816: }
817: PetscViewerASCIIPushTab(viewer);
818: KSPView(mglevels[i]->smoothd,viewer);
819: PetscViewerASCIIPopTab(viewer);
820: if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
821: PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");
822: } else if (i) {
823: PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);
824: PetscViewerASCIIPushTab(viewer);
825: KSPView(mglevels[i]->smoothu,viewer);
826: PetscViewerASCIIPopTab(viewer);
827: }
828: if (i && mglevels[i]->cr) {
829: PetscViewerASCIIPrintf(viewer,"CR solver on level %D -------------------------------\n",i);
830: PetscViewerASCIIPushTab(viewer);
831: KSPView(mglevels[i]->cr,viewer);
832: PetscViewerASCIIPopTab(viewer);
833: }
834: }
835: } else if (isbinary) {
836: for (i=levels-1; i>=0; i--) {
837: KSPView(mglevels[i]->smoothd,viewer);
838: if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) {
839: KSPView(mglevels[i]->smoothu,viewer);
840: }
841: }
842: } else if (isdraw) {
843: PetscDraw draw;
844: PetscReal x,w,y,bottom,th;
845: PetscViewerDrawGetDraw(viewer,0,&draw);
846: PetscDrawGetCurrentPoint(draw,&x,&y);
847: PetscDrawStringGetSize(draw,NULL,&th);
848: bottom = y - th;
849: for (i=levels-1; i>=0; i--) {
850: if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
851: PetscDrawPushCurrentPoint(draw,x,bottom);
852: KSPView(mglevels[i]->smoothd,viewer);
853: PetscDrawPopCurrentPoint(draw);
854: } else {
855: w = 0.5*PetscMin(1.0-x,x);
856: PetscDrawPushCurrentPoint(draw,x+w,bottom);
857: KSPView(mglevels[i]->smoothd,viewer);
858: PetscDrawPopCurrentPoint(draw);
859: PetscDrawPushCurrentPoint(draw,x-w,bottom);
860: KSPView(mglevels[i]->smoothu,viewer);
861: PetscDrawPopCurrentPoint(draw);
862: }
863: PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);
864: bottom -= th;
865: }
866: }
867: return(0);
868: }
870: #include <petsc/private/kspimpl.h>
872: /*
873: Calls setup for the KSP on each level
874: */
875: PetscErrorCode PCSetUp_MG(PC pc)
876: {
877: PC_MG *mg = (PC_MG*)pc->data;
878: PC_MG_Levels **mglevels = mg->levels;
880: PetscInt i,n;
881: PC cpc;
882: PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE;
883: Mat dA,dB;
884: Vec tvec;
885: DM *dms;
886: PetscViewer viewer = NULL;
887: PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE;
890: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up");
891: n = mglevels[0]->levels;
892: /* FIX: Move this to PCSetFromOptions_MG? */
893: if (mg->usedmfornumberoflevels) {
894: PetscInt levels;
895: DMGetRefineLevel(pc->dm,&levels);
896: levels++;
897: if (levels > n) { /* the problem is now being solved on a finer grid */
898: PCMGSetLevels(pc,levels,NULL);
899: n = levels;
900: PCSetFromOptions(pc); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
901: mglevels = mg->levels;
902: }
903: }
904: KSPGetPC(mglevels[0]->smoothd,&cpc);
907: /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
908: /* so use those from global PC */
909: /* Is this what we always want? What if user wants to keep old one? */
910: KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);
911: if (opsset) {
912: Mat mmat;
913: KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);
914: if (mmat == pc->pmat) opsset = PETSC_FALSE;
915: }
917: /* Create CR solvers */
918: PCMGGetAdaptCR(pc, &doCR);
919: if (doCR) {
920: const char *prefix;
922: PCGetOptionsPrefix(pc, &prefix);
923: for (i = 1; i < n; ++i) {
924: PC ipc, cr;
925: char crprefix[128];
927: KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr);
928: KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE);
929: PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i);
930: KSPSetOptionsPrefix(mglevels[i]->cr, prefix);
931: PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level);
932: KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV);
933: KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL);
934: KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED);
935: KSPGetPC(mglevels[i]->cr, &ipc);
937: PCSetType(ipc, PCCOMPOSITE);
938: PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE);
939: PCCompositeAddPCType(ipc, PCSOR);
940: CreateCR_Private(pc, i, &cr);
941: PCCompositeAddPC(ipc, cr);
942: PCDestroy(&cr);
944: KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd);
945: KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE);
946: PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i);
947: KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix);
948: PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr);
949: }
950: }
952: if (!opsset) {
953: PCGetUseAmat(pc,&use_amat);
954: if (use_amat) {
955: PetscInfo(pc,"Using outer operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
956: KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);
957: } else {
958: PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");
959: KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);
960: }
961: }
963: for (i=n-1; i>0; i--) {
964: if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
965: missinginterpolate = PETSC_TRUE;
966: continue;
967: }
968: }
970: KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);
971: if (dA == dB) dAeqdB = PETSC_TRUE;
972: if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) {
973: needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */
974: }
977: /*
978: Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
979: Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
980: */
981: if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
982: /* construct the interpolation from the DMs */
983: Mat p;
984: Vec rscale;
985: PetscMalloc1(n,&dms);
986: dms[n-1] = pc->dm;
987: /* Separately create them so we do not get DMKSP interference between levels */
988: for (i=n-2; i>-1; i--) {DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);}
989: /*
990: Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level.
991: Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options.
992: But it is safe to use -dm_mat_type.
994: The mat type should not be hardcoded like this, we need to find a better way.
995: DMSetMatType(dms[0],MATAIJ);
996: */
997: for (i=n-2; i>-1; i--) {
998: DMKSP kdm;
999: PetscBool dmhasrestrict, dmhasinject;
1000: KSPSetDM(mglevels[i]->smoothd,dms[i]);
1001: if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);}
1002: if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1003: KSPSetDM(mglevels[i]->smoothu,dms[i]);
1004: if (!needRestricts) {KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);}
1005: }
1006: if (mglevels[i]->cr) {
1007: KSPSetDM(mglevels[i]->cr,dms[i]);
1008: if (!needRestricts) {KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE);}
1009: }
1010: DMGetDMKSPWrite(dms[i],&kdm);
1011: /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take
1012: * a bitwise OR of computing the matrix, RHS, and initial iterate. */
1013: kdm->ops->computerhs = NULL;
1014: kdm->rhsctx = NULL;
1015: if (!mglevels[i+1]->interpolate) {
1016: DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);
1017: PCMGSetInterpolation(pc,i+1,p);
1018: if (rscale) {PCMGSetRScale(pc,i+1,rscale);}
1019: VecDestroy(&rscale);
1020: MatDestroy(&p);
1021: }
1022: DMHasCreateRestriction(dms[i],&dmhasrestrict);
1023: if (dmhasrestrict && !mglevels[i+1]->restrct){
1024: DMCreateRestriction(dms[i],dms[i+1],&p);
1025: PCMGSetRestriction(pc,i+1,p);
1026: MatDestroy(&p);
1027: }
1028: DMHasCreateInjection(dms[i],&dmhasinject);
1029: if (dmhasinject && !mglevels[i+1]->inject){
1030: DMCreateInjection(dms[i],dms[i+1],&p);
1031: PCMGSetInjection(pc,i+1,p);
1032: MatDestroy(&p);
1033: }
1034: }
1036: for (i=n-2; i>-1; i--) {DMDestroy(&dms[i]);}
1037: PetscFree(dms);
1038: }
1040: if (pc->dm && !pc->setupcalled) {
1041: /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
1042: KSPSetDM(mglevels[n-1]->smoothd,pc->dm);
1043: KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);
1044: if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) {
1045: KSPSetDM(mglevels[n-1]->smoothu,pc->dm);
1046: KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);
1047: }
1048: if (mglevels[n-1]->cr) {
1049: KSPSetDM(mglevels[n-1]->cr,pc->dm);
1050: KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE);
1051: }
1052: }
1054: if (mg->galerkin < PC_MG_GALERKIN_NONE) {
1055: Mat A,B;
1056: PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE;
1057: MatReuse reuse = MAT_INITIAL_MATRIX;
1059: if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE;
1060: if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE;
1061: if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
1062: for (i=n-2; i>-1; i--) {
1063: 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");
1064: if (!mglevels[i+1]->interpolate) {
1065: PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);
1066: }
1067: if (!mglevels[i+1]->restrct) {
1068: PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);
1069: }
1070: if (reuse == MAT_REUSE_MATRIX) {
1071: KSPGetOperators(mglevels[i]->smoothd,&A,&B);
1072: }
1073: if (doA) {
1074: MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);
1075: }
1076: if (doB) {
1077: MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);
1078: }
1079: /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
1080: if (!doA && dAeqdB) {
1081: if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)B);}
1082: A = B;
1083: } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
1084: KSPGetOperators(mglevels[i]->smoothd,&A,NULL);
1085: PetscObjectReference((PetscObject)A);
1086: }
1087: if (!doB && dAeqdB) {
1088: if (reuse == MAT_INITIAL_MATRIX) {PetscObjectReference((PetscObject)A);}
1089: B = A;
1090: } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
1091: KSPGetOperators(mglevels[i]->smoothd,NULL,&B);
1092: PetscObjectReference((PetscObject)B);
1093: }
1094: if (reuse == MAT_INITIAL_MATRIX) {
1095: KSPSetOperators(mglevels[i]->smoothd,A,B);
1096: PetscObjectDereference((PetscObject)A);
1097: PetscObjectDereference((PetscObject)B);
1098: }
1099: dA = A;
1100: dB = B;
1101: }
1102: }
1105: /* Adapt interpolation matrices */
1106: if (mg->adaptInterpolation) {
1107: mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */
1108: for (i = 0; i < n; ++i) {
1109: PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);
1110: if (i) {PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mg->Nc, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace);}
1111: }
1112: for (i = n-2; i > -1; --i) {
1113: PCMGRecomputeLevelOperators_Internal(pc, i);
1114: }
1115: }
1117: if (needRestricts && pc->dm) {
1118: for (i=n-2; i>=0; i--) {
1119: DM dmfine,dmcoarse;
1120: Mat Restrict,Inject;
1121: Vec rscale;
1122: KSPGetDM(mglevels[i+1]->smoothd,&dmfine);
1123: KSPGetDM(mglevels[i]->smoothd,&dmcoarse);
1124: PCMGGetRestriction(pc,i+1,&Restrict);
1125: PCMGGetRScale(pc,i+1,&rscale);
1126: PCMGGetInjection(pc,i+1,&Inject);
1127: DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);
1128: }
1129: }
1131: if (!pc->setupcalled) {
1132: for (i=0; i<n; i++) {
1133: KSPSetFromOptions(mglevels[i]->smoothd);
1134: }
1135: for (i=1; i<n; i++) {
1136: if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) {
1137: KSPSetFromOptions(mglevels[i]->smoothu);
1138: }
1139: if (mglevels[i]->cr) {
1140: KSPSetFromOptions(mglevels[i]->cr);
1141: }
1142: }
1143: /* insure that if either interpolation or restriction is set the other other one is set */
1144: for (i=1; i<n; i++) {
1145: PCMGGetInterpolation(pc,i,NULL);
1146: PCMGGetRestriction(pc,i,NULL);
1147: }
1148: for (i=0; i<n-1; i++) {
1149: if (!mglevels[i]->b) {
1150: Vec *vec;
1151: KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);
1152: PCMGSetRhs(pc,i,*vec);
1153: VecDestroy(vec);
1154: PetscFree(vec);
1155: }
1156: if (!mglevels[i]->r && i) {
1157: VecDuplicate(mglevels[i]->b,&tvec);
1158: PCMGSetR(pc,i,tvec);
1159: VecDestroy(&tvec);
1160: }
1161: if (!mglevels[i]->x) {
1162: VecDuplicate(mglevels[i]->b,&tvec);
1163: PCMGSetX(pc,i,tvec);
1164: VecDestroy(&tvec);
1165: }
1166: if (doCR) {
1167: VecDuplicate(mglevels[i]->b,&mglevels[i]->crx);
1168: VecDuplicate(mglevels[i]->b,&mglevels[i]->crb);
1169: VecZeroEntries(mglevels[i]->crb);
1170: }
1171: }
1172: if (n != 1 && !mglevels[n-1]->r) {
1173: /* PCMGSetR() on the finest level if user did not supply it */
1174: Vec *vec;
1175: KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);
1176: PCMGSetR(pc,n-1,*vec);
1177: VecDestroy(vec);
1178: PetscFree(vec);
1179: }
1180: if (doCR) {
1181: VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx);
1182: VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb);
1183: VecZeroEntries(mglevels[n-1]->crb);
1184: }
1185: }
1187: if (pc->dm) {
1188: /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
1189: for (i=0; i<n-1; i++) {
1190: if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
1191: }
1192: }
1194: for (i=1; i<n; i++) {
1195: if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){
1196: /* if doing only down then initial guess is zero */
1197: KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);
1198: }
1199: if (mglevels[i]->cr) {KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);}
1200: if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1201: KSPSetUp(mglevels[i]->smoothd);
1202: if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1203: pc->failedreason = PC_SUBPC_ERROR;
1204: }
1205: if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1206: if (!mglevels[i]->residual) {
1207: Mat mat;
1208: KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);
1209: PCMGSetResidual(pc,i,PCMGResidualDefault,mat);
1210: }
1211: if (!mglevels[i]->residualtranspose) {
1212: Mat mat;
1213: KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);
1214: PCMGSetResidualTranspose(pc,i,PCMGResidualTransposeDefault,mat);
1215: }
1216: }
1217: for (i=1; i<n; i++) {
1218: if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
1219: Mat downmat,downpmat;
1221: /* check if operators have been set for up, if not use down operators to set them */
1222: KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);
1223: if (!opsset) {
1224: KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
1225: KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);
1226: }
1228: KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);
1229: if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1230: KSPSetUp(mglevels[i]->smoothu);
1231: if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) {
1232: pc->failedreason = PC_SUBPC_ERROR;
1233: }
1234: if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1235: }
1236: if (mglevels[i]->cr) {
1237: Mat downmat,downpmat;
1239: /* check if operators have been set for up, if not use down operators to set them */
1240: KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL);
1241: if (!opsset) {
1242: KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);
1243: KSPSetOperators(mglevels[i]->cr,downmat,downpmat);
1244: }
1246: KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);
1247: if (mglevels[i]->eventsmoothsetup) {PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1248: KSPSetUp(mglevels[i]->cr);
1249: if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) {
1250: pc->failedreason = PC_SUBPC_ERROR;
1251: }
1252: if (mglevels[i]->eventsmoothsetup) {PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);}
1253: }
1254: }
1256: if (mglevels[0]->eventsmoothsetup) {PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);}
1257: KSPSetUp(mglevels[0]->smoothd);
1258: if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) {
1259: pc->failedreason = PC_SUBPC_ERROR;
1260: }
1261: if (mglevels[0]->eventsmoothsetup) {PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);}
1263: /*
1264: Dump the interpolation/restriction matrices plus the
1265: Jacobian/stiffness on each level. This allows MATLAB users to
1266: easily check if the Galerkin condition A_c = R A_f R^T is satisfied.
1268: Only support one or the other at the same time.
1269: */
1270: #if defined(PETSC_USE_SOCKET_VIEWER)
1271: PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);
1272: if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
1273: dump = PETSC_FALSE;
1274: #endif
1275: PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);
1276: if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));
1278: if (viewer) {
1279: for (i=1; i<n; i++) {
1280: MatView(mglevels[i]->restrct,viewer);
1281: }
1282: for (i=0; i<n; i++) {
1283: KSPGetPC(mglevels[i]->smoothd,&pc);
1284: MatView(pc->mat,viewer);
1285: }
1286: }
1287: return(0);
1288: }
1290: /* -------------------------------------------------------------------------------------*/
1292: PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
1293: {
1294: PC_MG *mg = (PC_MG *) pc->data;
1297: *levels = mg->nlevels;
1298: return(0);
1299: }
1301: /*@
1302: PCMGGetLevels - Gets the number of levels to use with MG.
1304: Not Collective
1306: Input Parameter:
1307: . pc - the preconditioner context
1309: Output parameter:
1310: . levels - the number of levels
1312: Level: advanced
1314: .seealso: PCMGSetLevels()
1315: @*/
1316: PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels)
1317: {
1323: *levels = 0;
1324: PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));
1325: return(0);
1326: }
1328: /*@
1329: PCMGSetType - Determines the form of multigrid to use:
1330: multiplicative, additive, full, or the Kaskade algorithm.
1332: Logically Collective on PC
1334: Input Parameters:
1335: + pc - the preconditioner context
1336: - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,
1337: PC_MG_FULL, PC_MG_KASKADE
1339: Options Database Key:
1340: . -pc_mg_type <form> - Sets <form>, one of multiplicative,
1341: additive, full, kaskade
1343: Level: advanced
1345: .seealso: PCMGSetLevels()
1346: @*/
1347: PetscErrorCode PCMGSetType(PC pc,PCMGType form)
1348: {
1349: PC_MG *mg = (PC_MG*)pc->data;
1354: mg->am = form;
1355: if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1356: else pc->ops->applyrichardson = NULL;
1357: return(0);
1358: }
1360: /*@
1361: PCMGGetType - Determines the form of multigrid to use:
1362: multiplicative, additive, full, or the Kaskade algorithm.
1364: Logically Collective on PC
1366: Input Parameter:
1367: . pc - the preconditioner context
1369: Output Parameter:
1370: . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE
1373: Level: advanced
1375: .seealso: PCMGSetLevels()
1376: @*/
1377: PetscErrorCode PCMGGetType(PC pc,PCMGType *type)
1378: {
1379: PC_MG *mg = (PC_MG*)pc->data;
1383: *type = mg->am;
1384: return(0);
1385: }
1387: /*@
1388: PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more
1389: complicated cycling.
1391: Logically Collective on PC
1393: Input Parameters:
1394: + pc - the multigrid context
1395: - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W
1397: Options Database Key:
1398: . -pc_mg_cycle_type <v,w> - provide the cycle desired
1400: Level: advanced
1402: .seealso: PCMGSetCycleTypeOnLevel()
1403: @*/
1404: PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n)
1405: {
1406: PC_MG *mg = (PC_MG*)pc->data;
1407: PC_MG_Levels **mglevels = mg->levels;
1408: PetscInt i,levels;
1413: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1414: levels = mglevels[0]->levels;
1415: for (i=0; i<levels; i++) mglevels[i]->cycles = n;
1416: return(0);
1417: }
1419: /*@
1420: PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
1421: of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used
1423: Logically Collective on PC
1425: Input Parameters:
1426: + pc - the multigrid context
1427: - n - number of cycles (default is 1)
1429: Options Database Key:
1430: . -pc_mg_multiplicative_cycles n
1432: Level: advanced
1434: Notes:
1435: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType()
1437: .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType()
1438: @*/
1439: PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n)
1440: {
1441: PC_MG *mg = (PC_MG*)pc->data;
1446: mg->cyclesperpcapply = n;
1447: return(0);
1448: }
1450: PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use)
1451: {
1452: PC_MG *mg = (PC_MG*)pc->data;
1455: mg->galerkin = use;
1456: return(0);
1457: }
1459: /*@
1460: PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
1461: finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i
1463: Logically Collective on PC
1465: Input Parameters:
1466: + pc - the multigrid context
1467: - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE
1469: Options Database Key:
1470: . -pc_mg_galerkin <both,pmat,mat,none>
1472: Level: intermediate
1474: Notes:
1475: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not
1476: use the PCMG construction of the coarser grids.
1478: .seealso: PCMGGetGalerkin(), PCMGGalerkinType
1480: @*/
1481: PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use)
1482: {
1487: PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));
1488: return(0);
1489: }
1491: /*@
1492: PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e.
1493: A_i-1 = r_i * A_i * p_i
1495: Not Collective
1497: Input Parameter:
1498: . pc - the multigrid context
1500: Output Parameter:
1501: . galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL
1503: Level: intermediate
1505: .seealso: PCMGSetGalerkin(), PCMGGalerkinType
1507: @*/
1508: PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin)
1509: {
1510: PC_MG *mg = (PC_MG*)pc->data;
1514: *galerkin = mg->galerkin;
1515: return(0);
1516: }
1518: PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1519: {
1520: PC_MG *mg = (PC_MG *) pc->data;
1523: mg->adaptInterpolation = adapt;
1524: return(0);
1525: }
1527: PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1528: {
1529: PC_MG *mg = (PC_MG *) pc->data;
1532: *adapt = mg->adaptInterpolation;
1533: return(0);
1534: }
1536: PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr)
1537: {
1538: PC_MG *mg = (PC_MG *) pc->data;
1541: mg->compatibleRelaxation = cr;
1542: return(0);
1543: }
1545: PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr)
1546: {
1547: PC_MG *mg = (PC_MG *) pc->data;
1550: *cr = mg->compatibleRelaxation;
1551: return(0);
1552: }
1554: /*@
1555: PCMGSetAdaptInterpolation - Adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated.
1557: Logically Collective on PC
1559: Input Parameters:
1560: + pc - the multigrid context
1561: - adapt - flag for adaptation of the interpolator
1563: Options Database Keys:
1564: + -pc_mg_adapt_interp - Turn on adaptation
1565: . -pc_mg_adapt_interp_n <int> - The number of modes to use, should be divisible by dimension
1566: - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector
1568: Level: intermediate
1570: .keywords: MG, set, Galerkin
1571: .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1572: @*/
1573: PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1574: {
1579: PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));
1580: return(0);
1581: }
1583: /*@
1584: 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.
1586: Not collective
1588: Input Parameter:
1589: . pc - the multigrid context
1591: Output Parameter:
1592: . adapt - flag for adaptation of the interpolator
1594: Level: intermediate
1596: .keywords: MG, set, Galerkin
1597: .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1598: @*/
1599: PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1600: {
1606: PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));
1607: return(0);
1608: }
1610: /*@
1611: PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation.
1613: Logically Collective on PC
1615: Input Parameters:
1616: + pc - the multigrid context
1617: - cr - flag for compatible relaxation
1619: Options Database Keys:
1620: . -pc_mg_adapt_cr - Turn on compatible relaxation
1622: Level: intermediate
1624: .keywords: MG, set, Galerkin
1625: .seealso: PCMGGetAdaptCR(), PCMGSetAdaptInterpolation(), PCMGSetGalerkin()
1626: @*/
1627: PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr)
1628: {
1633: PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr));
1634: return(0);
1635: }
1637: /*@
1638: PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation.
1640: Not collective
1642: Input Parameter:
1643: . pc - the multigrid context
1645: Output Parameter:
1646: . cr - flag for compatible relaxaion
1648: Level: intermediate
1650: .keywords: MG, set, Galerkin
1651: .seealso: PCMGSetAdaptCR(), PCMGGetAdaptInterpolation(), PCMGSetGalerkin()
1652: @*/
1653: PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr)
1654: {
1660: PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr));
1661: return(0);
1662: }
1664: /*@
1665: PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1666: on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of
1667: pre- and post-smoothing steps.
1669: Logically Collective on PC
1671: Input Parameters:
1672: + mg - the multigrid context
1673: - n - the number of smoothing steps
1675: Options Database Key:
1676: . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps
1678: Level: advanced
1680: Notes:
1681: this does not set a value on the coarsest grid, since we assume that
1682: there is no separate smooth up on the coarsest grid.
1684: .seealso: PCMGSetDistinctSmoothUp()
1685: @*/
1686: PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n)
1687: {
1688: PC_MG *mg = (PC_MG*)pc->data;
1689: PC_MG_Levels **mglevels = mg->levels;
1691: PetscInt i,levels;
1696: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1697: levels = mglevels[0]->levels;
1699: for (i=1; i<levels; i++) {
1700: KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1701: KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);
1702: mg->default_smoothu = n;
1703: mg->default_smoothd = n;
1704: }
1705: return(0);
1706: }
1708: /*@
1709: PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels
1710: and adds the suffix _up to the options name
1712: Logically Collective on PC
1714: Input Parameters:
1715: . pc - the preconditioner context
1717: Options Database Key:
1718: . -pc_mg_distinct_smoothup
1720: Level: advanced
1722: Notes:
1723: this does not set a value on the coarsest grid, since we assume that
1724: there is no separate smooth up on the coarsest grid.
1726: .seealso: PCMGSetNumberSmooth()
1727: @*/
1728: PetscErrorCode PCMGSetDistinctSmoothUp(PC pc)
1729: {
1730: PC_MG *mg = (PC_MG*)pc->data;
1731: PC_MG_Levels **mglevels = mg->levels;
1733: PetscInt i,levels;
1734: KSP subksp;
1738: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling");
1739: levels = mglevels[0]->levels;
1741: for (i=1; i<levels; i++) {
1742: const char *prefix = NULL;
1743: /* make sure smoother up and down are different */
1744: PCMGGetSmootherUp(pc,i,&subksp);
1745: KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);
1746: KSPSetOptionsPrefix(subksp,prefix);
1747: KSPAppendOptionsPrefix(subksp,"up_");
1748: }
1749: return(0);
1750: }
1752: /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1753: PetscErrorCode PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[])
1754: {
1755: PC_MG *mg = (PC_MG*)pc->data;
1756: PC_MG_Levels **mglevels = mg->levels;
1757: Mat *mat;
1758: PetscInt l;
1762: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1763: PetscMalloc1(mg->nlevels,&mat);
1764: for (l=1; l< mg->nlevels; l++) {
1765: mat[l-1] = mglevels[l]->interpolate;
1766: PetscObjectReference((PetscObject)mat[l-1]);
1767: }
1768: *num_levels = mg->nlevels;
1769: *interpolations = mat;
1770: return(0);
1771: }
1773: /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1774: PetscErrorCode PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[])
1775: {
1776: PC_MG *mg = (PC_MG*)pc->data;
1777: PC_MG_Levels **mglevels = mg->levels;
1778: PetscInt l;
1779: Mat *mat;
1783: if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling");
1784: PetscMalloc1(mg->nlevels,&mat);
1785: for (l=0; l<mg->nlevels-1; l++) {
1786: KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));
1787: PetscObjectReference((PetscObject)mat[l]);
1788: }
1789: *num_levels = mg->nlevels;
1790: *coarseOperators = mat;
1791: return(0);
1792: }
1794: /*@C
1795: PCMGRegisterCoarseSpaceConstructor - Adds a method to the PCMG package for coarse space construction.
1797: Not collective
1799: Input Parameters:
1800: + name - name of the constructor
1801: - function - constructor routine
1803: Notes:
1804: Calling sequence for the routine:
1805: $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1806: $ pc - The PC object
1807: $ l - The multigrid level, 0 is the coarse level
1808: $ dm - The DM for this level
1809: $ smooth - The level smoother
1810: $ Nc - The size of the coarse space
1811: $ initGuess - Basis for an initial guess for the space
1812: $ coarseSp - A basis for the computed coarse space
1814: Level: advanced
1816: .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister()
1817: @*/
1818: PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1819: {
1823: PCInitializePackage();
1824: PetscFunctionListAdd(&PCMGCoarseList,name,function);
1825: return(0);
1826: }
1828: /*@C
1829: PCMGGetCoarseSpaceConstructor - Returns the given coarse space construction method.
1831: Not collective
1833: Input Parameter:
1834: . name - name of the constructor
1836: Output Parameter:
1837: . function - constructor routine
1839: Notes:
1840: Calling sequence for the routine:
1841: $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp)
1842: $ pc - The PC object
1843: $ l - The multigrid level, 0 is the coarse level
1844: $ dm - The DM for this level
1845: $ smooth - The level smoother
1846: $ Nc - The size of the coarse space
1847: $ initGuess - Basis for an initial guess for the space
1848: $ coarseSp - A basis for the computed coarse space
1850: Level: advanced
1852: .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister()
1853: @*/
1854: PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **))
1855: {
1859: PetscFunctionListFind(PCMGCoarseList,name,function);
1860: return(0);
1861: }
1863: /* ----------------------------------------------------------------------------------------*/
1865: /*MC
1866: PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional
1867: information about the coarser grid matrices and restriction/interpolation operators.
1869: Options Database Keys:
1870: + -pc_mg_levels <nlevels> - number of levels including finest
1871: . -pc_mg_cycle_type <v,w> - provide the cycle desired
1872: . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default
1873: . -pc_mg_log - log information about time spent on each level of the solver
1874: . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
1875: . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R'
1876: . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1)
1877: . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices
1878: to the Socket viewer for reading from MATLAB.
1879: - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices
1880: to the binary output file called binaryoutput
1882: Notes:
1883: 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
1885: When run with a single level the smoother options are used on that level NOT the coarse grid solver options
1887: When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
1888: is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
1889: (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
1890: residual is computed at the end of each cycle.
1892: Level: intermediate
1894: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE
1895: PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(),
1896: PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
1897: PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
1898: PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
1899: M*/
1901: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1902: {
1903: PC_MG *mg;
1907: PetscNewLog(pc,&mg);
1908: pc->data = (void*)mg;
1909: mg->nlevels = -1;
1910: mg->am = PC_MG_MULTIPLICATIVE;
1911: mg->galerkin = PC_MG_GALERKIN_NONE;
1912: mg->adaptInterpolation = PETSC_FALSE;
1913: mg->Nc = -1;
1914: mg->eigenvalue = -1;
1916: pc->useAmat = PETSC_TRUE;
1918: pc->ops->apply = PCApply_MG;
1919: pc->ops->applytranspose = PCApplyTranspose_MG;
1920: pc->ops->matapply = PCMatApply_MG;
1921: pc->ops->setup = PCSetUp_MG;
1922: pc->ops->reset = PCReset_MG;
1923: pc->ops->destroy = PCDestroy_MG;
1924: pc->ops->setfromoptions = PCSetFromOptions_MG;
1925: pc->ops->view = PCView_MG;
1927: PetscObjectComposedDataRegister(&mg->eigenvalue);
1928: PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);
1929: PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);
1930: PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);
1931: PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);
1932: PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);
1933: PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);
1934: PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);
1935: PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG);
1936: PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG);
1937: return(0);
1938: }