Actual source code: mg.c

  1: /*
  2:     Defines the multigrid preconditioner interface.
  3: */
  4: #include <petsc/private/pcmgimpl.h>
  5: #include <petsc/private/kspimpl.h>
  6: #include <petscdm.h>
  7: PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC, PetscBool *);

  9: /*
 10:    Contains the list of registered coarse space construction routines
 11: */
 12: PetscFunctionList PCMGCoarseList = NULL;

 14: PetscErrorCode PCMGMCycle_Private(PC pc, PC_MG_Levels **mglevelsin, PetscBool transpose, PetscBool matapp, PCRichardsonConvergedReason *reason)
 15: {
 16:   PC_MG        *mg = (PC_MG *)pc->data;
 17:   PC_MG_Levels *mgc, *mglevels = *mglevelsin;
 18:   PetscInt      cycles = (mglevels->level == 1) ? 1 : mglevels->cycles;

 20:   PetscFunctionBegin;
 21:   if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve, 0, 0, 0, 0));
 22:   if (!transpose) {
 23:     if (matapp) {
 24:       PetscCall(KSPMatSolve(mglevels->smoothd, mglevels->B, mglevels->X)); /* pre-smooth */
 25:       PetscCall(KSPCheckSolve(mglevels->smoothd, pc, NULL));
 26:     } else {
 27:       PetscCall(KSPSolve(mglevels->smoothd, mglevels->b, mglevels->x)); /* pre-smooth */
 28:       PetscCall(KSPCheckSolve(mglevels->smoothd, pc, mglevels->x));
 29:     }
 30:   } else {
 31:     PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported");
 32:     PetscCall(KSPSolveTranspose(mglevels->smoothu, mglevels->b, mglevels->x)); /* transpose of post-smooth */
 33:     PetscCall(KSPCheckSolve(mglevels->smoothu, pc, mglevels->x));
 34:   }
 35:   if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve, 0, 0, 0, 0));
 36:   if (mglevels->level) { /* not the coarsest grid */
 37:     if (mglevels->eventresidual) PetscCall(PetscLogEventBegin(mglevels->eventresidual, 0, 0, 0, 0));
 38:     if (matapp && !mglevels->R) PetscCall(MatDuplicate(mglevels->B, MAT_DO_NOT_COPY_VALUES, &mglevels->R));
 39:     if (!transpose) {
 40:       if (matapp) PetscCall((*mglevels->matresidual)(mglevels->A, mglevels->B, mglevels->X, mglevels->R));
 41:       else PetscCall((*mglevels->residual)(mglevels->A, mglevels->b, mglevels->x, mglevels->r));
 42:     } else {
 43:       if (matapp) PetscCall((*mglevels->matresidualtranspose)(mglevels->A, mglevels->B, mglevels->X, mglevels->R));
 44:       else PetscCall((*mglevels->residualtranspose)(mglevels->A, mglevels->b, mglevels->x, mglevels->r));
 45:     }
 46:     if (mglevels->eventresidual) PetscCall(PetscLogEventEnd(mglevels->eventresidual, 0, 0, 0, 0));

 48:     /* if on finest level and have convergence criteria set */
 49:     if (mglevels->level == mglevels->levels - 1 && mg->ttol && reason) {
 50:       PetscReal rnorm;

 52:       PetscCall(VecNorm(mglevels->r, NORM_2, &rnorm));
 53:       if (rnorm <= mg->ttol) {
 54:         if (rnorm < mg->abstol) {
 55:           *reason = PCRICHARDSON_CONVERGED_ATOL;
 56:           PetscCall(PetscInfo(pc, "Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n", (double)rnorm, (double)mg->abstol));
 57:         } else {
 58:           *reason = PCRICHARDSON_CONVERGED_RTOL;
 59:           PetscCall(PetscInfo(pc, "Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n", (double)rnorm, (double)mg->ttol));
 60:         }
 61:         PetscFunctionReturn(PETSC_SUCCESS);
 62:       }
 63:     }

 65:     mgc = *(mglevelsin - 1);
 66:     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict, 0, 0, 0, 0));
 67:     if (!transpose) {
 68:       if (matapp) PetscCall(MatMatRestrict(mglevels->restrct, mglevels->R, &mgc->B));
 69:       else PetscCall(MatRestrict(mglevels->restrct, mglevels->r, mgc->b));
 70:     } else {
 71:       if (matapp) PetscCall(MatMatRestrict(mglevels->interpolate, mglevels->R, &mgc->B));
 72:       else PetscCall(MatRestrict(mglevels->interpolate, mglevels->r, mgc->b));
 73:     }
 74:     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict, 0, 0, 0, 0));
 75:     if (matapp) {
 76:       if (!mgc->X) PetscCall(MatDuplicate(mgc->B, MAT_DO_NOT_COPY_VALUES, &mgc->X));
 77:       else {
 78:         PetscCall(MatZeroEntries(mgc->X));
 79:       }
 80:     } else {
 81:       PetscCall(VecZeroEntries(mgc->x));
 82:     }
 83:     while (cycles--) PetscCall(PCMGMCycle_Private(pc, mglevelsin - 1, transpose, matapp, reason));
 84:     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict, 0, 0, 0, 0));
 85:     if (!transpose) {
 86:       if (matapp) PetscCall(MatMatInterpolateAdd(mglevels->interpolate, mgc->X, mglevels->X, &mglevels->X));
 87:       else PetscCall(MatInterpolateAdd(mglevels->interpolate, mgc->x, mglevels->x, mglevels->x));
 88:     } else {
 89:       PetscCall(MatInterpolateAdd(mglevels->restrct, mgc->x, mglevels->x, mglevels->x));
 90:     }
 91:     if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict, 0, 0, 0, 0));
 92:     if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve, 0, 0, 0, 0));
 93:     if (!transpose) {
 94:       if (matapp) {
 95:         PetscCall(KSPMatSolve(mglevels->smoothu, mglevels->B, mglevels->X)); /* post smooth */
 96:         PetscCall(KSPCheckSolve(mglevels->smoothu, pc, NULL));
 97:       } else {
 98:         PetscCall(KSPSolve(mglevels->smoothu, mglevels->b, mglevels->x)); /* post smooth */
 99:         PetscCall(KSPCheckSolve(mglevels->smoothu, pc, mglevels->x));
100:       }
101:     } else {
102:       PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported");
103:       PetscCall(KSPSolveTranspose(mglevels->smoothd, mglevels->b, mglevels->x)); /* post smooth */
104:       PetscCall(KSPCheckSolve(mglevels->smoothd, pc, mglevels->x));
105:     }
106:     if (mglevels->cr) {
107:       Mat crA;

109:       PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported");
110:       /* TODO Turn on copy and turn off noisy if we have an exact solution
111:       PetscCall(VecCopy(mglevels->x, mglevels->crx));
112:       PetscCall(VecCopy(mglevels->b, mglevels->crb)); */
113:       PetscCall(KSPGetOperators(mglevels->cr, &crA, NULL));
114:       PetscCall(KSPSetNoisy_Private(crA, mglevels->crx));
115:       PetscCall(KSPSolve(mglevels->cr, mglevels->crb, mglevels->crx)); /* compatible relaxation */
116:       PetscCall(KSPCheckSolve(mglevels->cr, pc, mglevels->crx));
117:     }
118:     if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve, 0, 0, 0, 0));
119:   }
120:   PetscFunctionReturn(PETSC_SUCCESS);
121: }

123: 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)
124: {
125:   PC_MG         *mg       = (PC_MG *)pc->data;
126:   PC_MG_Levels **mglevels = mg->levels;
127:   PC             tpc;
128:   PetscBool      changeu, changed;
129:   PetscInt       levels = mglevels[0]->levels, i;

131:   PetscFunctionBegin;
132:   /* When the DM is supplying the matrix then it will not exist until here */
133:   for (i = 0; i < levels; i++) {
134:     if (!mglevels[i]->A) {
135:       PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL));
136:       PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A));
137:     }
138:   }

140:   PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc));
141:   PetscCall(PCPreSolveChangeRHS(tpc, &changed));
142:   PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc));
143:   PetscCall(PCPreSolveChangeRHS(tpc, &changeu));
144:   if (!changed && !changeu) {
145:     PetscCall(VecDestroy(&mglevels[levels - 1]->b));
146:     mglevels[levels - 1]->b = b;
147:   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
148:     if (!mglevels[levels - 1]->b) {
149:       Vec *vec;

151:       PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL));
152:       mglevels[levels - 1]->b = *vec;
153:       PetscCall(PetscFree(vec));
154:     }
155:     PetscCall(VecCopy(b, mglevels[levels - 1]->b));
156:   }
157:   mglevels[levels - 1]->x = x;

159:   mg->rtol   = rtol;
160:   mg->abstol = abstol;
161:   mg->dtol   = dtol;
162:   if (rtol) {
163:     /* compute initial residual norm for relative convergence test */
164:     PetscReal rnorm;

166:     if (zeroguess) {
167:       PetscCall(VecNorm(b, NORM_2, &rnorm));
168:     } else {
169:       PetscCall((*mglevels[levels - 1]->residual)(mglevels[levels - 1]->A, b, x, w));
170:       PetscCall(VecNorm(w, NORM_2, &rnorm));
171:     }
172:     mg->ttol = PetscMax(rtol * rnorm, abstol);
173:   } else if (abstol) mg->ttol = abstol;
174:   else mg->ttol = 0.0;

176:   /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't
177:      stop prematurely due to small residual */
178:   for (i = 1; i < levels; i++) {
179:     PetscCall(KSPSetTolerances(mglevels[i]->smoothu, 0, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
180:     if (mglevels[i]->smoothu != mglevels[i]->smoothd) {
181:       /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */
182:       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE));
183:       PetscCall(KSPSetTolerances(mglevels[i]->smoothd, 0, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT));
184:     }
185:   }

187:   *reason = PCRICHARDSON_NOT_SET;
188:   for (i = 0; i < its; i++) {
189:     PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, PETSC_FALSE, PETSC_FALSE, reason));
190:     if (*reason) break;
191:   }
192:   if (*reason == PCRICHARDSON_NOT_SET) *reason = PCRICHARDSON_CONVERGED_ITS;
193:   *outits = i;
194:   if (!changed && !changeu) mglevels[levels - 1]->b = NULL;
195:   PetscFunctionReturn(PETSC_SUCCESS);
196: }

198: PetscErrorCode PCReset_MG(PC pc)
199: {
200:   PC_MG         *mg       = (PC_MG *)pc->data;
201:   PC_MG_Levels **mglevels = mg->levels;
202:   PetscInt       i, n;

204:   PetscFunctionBegin;
205:   if (mglevels) {
206:     n = mglevels[0]->levels;
207:     for (i = 0; i < n - 1; i++) {
208:       PetscCall(VecDestroy(&mglevels[i + 1]->r));
209:       PetscCall(VecDestroy(&mglevels[i]->b));
210:       PetscCall(VecDestroy(&mglevels[i]->x));
211:       PetscCall(MatDestroy(&mglevels[i + 1]->R));
212:       PetscCall(MatDestroy(&mglevels[i]->B));
213:       PetscCall(MatDestroy(&mglevels[i]->X));
214:       PetscCall(VecDestroy(&mglevels[i]->crx));
215:       PetscCall(VecDestroy(&mglevels[i]->crb));
216:       PetscCall(MatDestroy(&mglevels[i + 1]->restrct));
217:       PetscCall(MatDestroy(&mglevels[i + 1]->interpolate));
218:       PetscCall(MatDestroy(&mglevels[i + 1]->inject));
219:       PetscCall(VecDestroy(&mglevels[i + 1]->rscale));
220:     }
221:     PetscCall(VecDestroy(&mglevels[n - 1]->crx));
222:     PetscCall(VecDestroy(&mglevels[n - 1]->crb));
223:     /* this is not null only if the smoother on the finest level
224:        changes the rhs during PreSolve */
225:     PetscCall(VecDestroy(&mglevels[n - 1]->b));
226:     PetscCall(MatDestroy(&mglevels[n - 1]->B));

228:     for (i = 0; i < n; i++) {
229:       PetscCall(MatDestroy(&mglevels[i]->coarseSpace));
230:       PetscCall(MatDestroy(&mglevels[i]->A));
231:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPReset(mglevels[i]->smoothd));
232:       PetscCall(KSPReset(mglevels[i]->smoothu));
233:       if (mglevels[i]->cr) PetscCall(KSPReset(mglevels[i]->cr));
234:     }
235:     mg->Nc = 0;
236:   }
237:   PetscFunctionReturn(PETSC_SUCCESS);
238: }

240: /* Implementing CR

242: 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

244:   Inj^T (Inj Inj^T)^{-1} Inj

246: and if Inj is a VecScatter, as it is now in PETSc, we have

248:   Inj^T Inj

250: and

252:   S = I - Inj^T Inj

254: since

256:   Inj S = Inj - (Inj Inj^T) Inj = 0.

258: Brannick suggests

260:   A \to S^T A S  \qquad\mathrm{and}\qquad M \to S^T M S

262: 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

264:   M^{-1} A \to S M^{-1} A S

266: In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left.

268:   Check: || Inj P - I ||_F < tol
269:   Check: In general, Inj Inj^T = I
270: */

272: typedef struct {
273:   PC       mg;  /* The PCMG object */
274:   PetscInt l;   /* The multigrid level for this solver */
275:   Mat      Inj; /* The injection matrix */
276:   Mat      S;   /* I - Inj^T Inj */
277: } CRContext;

279: static PetscErrorCode CRSetup_Private(PC pc)
280: {
281:   CRContext *ctx;
282:   Mat        It;

284:   PetscFunctionBeginUser;
285:   PetscCall(PCShellGetContext(pc, &ctx));
286:   PetscCall(PCMGGetInjection(ctx->mg, ctx->l, &It));
287:   PetscCheck(It, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG");
288:   PetscCall(MatCreateTranspose(It, &ctx->Inj));
289:   PetscCall(MatCreateNormal(ctx->Inj, &ctx->S));
290:   PetscCall(MatScale(ctx->S, -1.0));
291:   PetscCall(MatShift(ctx->S, 1.0));
292:   PetscFunctionReturn(PETSC_SUCCESS);
293: }

295: static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y)
296: {
297:   CRContext *ctx;

299:   PetscFunctionBeginUser;
300:   PetscCall(PCShellGetContext(pc, &ctx));
301:   PetscCall(MatMult(ctx->S, x, y));
302:   PetscFunctionReturn(PETSC_SUCCESS);
303: }

305: static PetscErrorCode CRDestroy_Private(PC pc)
306: {
307:   CRContext *ctx;

309:   PetscFunctionBeginUser;
310:   PetscCall(PCShellGetContext(pc, &ctx));
311:   PetscCall(MatDestroy(&ctx->Inj));
312:   PetscCall(MatDestroy(&ctx->S));
313:   PetscCall(PetscFree(ctx));
314:   PetscCall(PCShellSetContext(pc, NULL));
315:   PetscFunctionReturn(PETSC_SUCCESS);
316: }

318: static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr)
319: {
320:   CRContext *ctx;

322:   PetscFunctionBeginUser;
323:   PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), cr));
324:   PetscCall(PetscObjectSetName((PetscObject)*cr, "S (complementary projector to injection)"));
325:   PetscCall(PetscCalloc1(1, &ctx));
326:   ctx->mg = pc;
327:   ctx->l  = l;
328:   PetscCall(PCSetType(*cr, PCSHELL));
329:   PetscCall(PCShellSetContext(*cr, ctx));
330:   PetscCall(PCShellSetApply(*cr, CRApply_Private));
331:   PetscCall(PCShellSetSetUp(*cr, CRSetup_Private));
332:   PetscCall(PCShellSetDestroy(*cr, CRDestroy_Private));
333:   PetscFunctionReturn(PETSC_SUCCESS);
334: }

336: PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions, const char[], const char[], const char *[], const char *[], PetscBool *);

338: PetscErrorCode PCMGSetLevels_MG(PC pc, PetscInt levels, MPI_Comm *comms)
339: {
340:   PC_MG         *mg = (PC_MG *)pc->data;
341:   MPI_Comm       comm;
342:   PC_MG_Levels **mglevels = mg->levels;
343:   PCMGType       mgtype   = mg->am;
344:   PetscInt       mgctype  = (PetscInt)PC_MG_CYCLE_V;
345:   PetscInt       i;
346:   PetscMPIInt    size;
347:   const char    *prefix;
348:   PC             ipc;
349:   PetscInt       n;

351:   PetscFunctionBegin;
354:   if (mg->nlevels == levels) PetscFunctionReturn(PETSC_SUCCESS);
355:   PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
356:   if (mglevels) {
357:     mgctype = mglevels[0]->cycles;
358:     /* changing the number of levels so free up the previous stuff */
359:     PetscCall(PCReset_MG(pc));
360:     n = mglevels[0]->levels;
361:     for (i = 0; i < n; i++) {
362:       if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPDestroy(&mglevels[i]->smoothd));
363:       PetscCall(KSPDestroy(&mglevels[i]->smoothu));
364:       PetscCall(KSPDestroy(&mglevels[i]->cr));
365:       PetscCall(PetscFree(mglevels[i]));
366:     }
367:     PetscCall(PetscFree(mg->levels));
368:   }

370:   mg->nlevels = levels;

372:   PetscCall(PetscMalloc1(levels, &mglevels));

374:   PetscCall(PCGetOptionsPrefix(pc, &prefix));

376:   mg->stageApply = 0;
377:   for (i = 0; i < levels; i++) {
378:     PetscCall(PetscNew(&mglevels[i]));

380:     mglevels[i]->level               = i;
381:     mglevels[i]->levels              = levels;
382:     mglevels[i]->cycles              = mgctype;
383:     mg->default_smoothu              = 2;
384:     mg->default_smoothd              = 2;
385:     mglevels[i]->eventsmoothsetup    = 0;
386:     mglevels[i]->eventsmoothsolve    = 0;
387:     mglevels[i]->eventresidual       = 0;
388:     mglevels[i]->eventinterprestrict = 0;

390:     if (comms) comm = comms[i];
391:     if (comm != MPI_COMM_NULL) {
392:       PetscCall(KSPCreate(comm, &mglevels[i]->smoothd));
393:       PetscCall(KSPSetNestLevel(mglevels[i]->smoothd, pc->kspnestlevel));
394:       PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->smoothd, pc->erroriffailure));
395:       PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd, (PetscObject)pc, levels - i));
396:       PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd, prefix));
397:       PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level));
398:       if (i == 0 && levels > 1) { // coarse grid
399:         PetscCall(KSPAppendOptionsPrefix(mglevels[0]->smoothd, "mg_coarse_"));

401:         /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */
402:         PetscCall(KSPSetType(mglevels[0]->smoothd, KSPPREONLY));
403:         PetscCall(KSPGetPC(mglevels[0]->smoothd, &ipc));
404:         PetscCallMPI(MPI_Comm_size(comm, &size));
405:         if (size > 1) PetscCall(PCSetType(ipc, PCREDUNDANT));
406:         else PetscCall(PCSetType(ipc, PCLU));
407:         PetscCall(PCFactorSetShiftType(ipc, MAT_SHIFT_INBLOCKS));
408:       } else {
409:         char tprefix[128];

411:         PetscCall(KSPSetType(mglevels[i]->smoothd, KSPCHEBYSHEV));
412:         PetscCall(KSPSetConvergenceTest(mglevels[i]->smoothd, KSPConvergedSkip, NULL, NULL));
413:         PetscCall(KSPSetNormType(mglevels[i]->smoothd, KSP_NORM_NONE));
414:         PetscCall(KSPGetPC(mglevels[i]->smoothd, &ipc));
415:         PetscCall(PCSetType(ipc, PCSOR));
416:         PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, mg->default_smoothd));

418:         if (i == levels - 1 && levels > 1) { // replace 'mg_finegrid_' with 'mg_levels_X_'
419:           PetscBool set;

421:           PetscCall(PetscOptionsFindPairPrefix_Private(((PetscObject)mglevels[i]->smoothd)->options, ((PetscObject)mglevels[i]->smoothd)->prefix, "-mg_fine_", NULL, NULL, &set));
422:           if (set) {
423:             if (prefix) PetscCall(PetscSNPrintf(tprefix, 128, "%smg_fine_", prefix));
424:             else PetscCall(PetscSNPrintf(tprefix, 128, "mg_fine_"));
425:             PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd, tprefix));
426:           } else {
427:             PetscCall(PetscSNPrintf(tprefix, 128, "mg_levels_%" PetscInt_FMT "_", i));
428:             PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd, tprefix));
429:           }
430:         } else {
431:           PetscCall(PetscSNPrintf(tprefix, 128, "mg_levels_%" PetscInt_FMT "_", i));
432:           PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd, tprefix));
433:         }
434:       }
435:     }
436:     mglevels[i]->smoothu = mglevels[i]->smoothd;
437:     mg->rtol             = 0.0;
438:     mg->abstol           = 0.0;
439:     mg->dtol             = 0.0;
440:     mg->ttol             = 0.0;
441:     mg->cyclesperpcapply = 1;
442:   }
443:   mg->levels = mglevels;
444:   PetscCall(PCMGSetType(pc, mgtype));
445:   PetscFunctionReturn(PETSC_SUCCESS);
446: }

448: /*@C
449:   PCMGSetLevels - Sets the number of levels to use with `PCMG`.
450:   Must be called before any other `PCMG` routine.

452:   Logically Collective

454:   Input Parameters:
455: + pc     - the preconditioner context
456: . levels - the number of levels
457: - comms  - optional communicators for each level; this is to allow solving the coarser problems
458:            on smaller sets of processes. For processes that are not included in the computation
459:            you must pass `MPI_COMM_NULL`. Use comms = `NULL` to specify that all processes
460:            should participate in each level of problem.

462:   Options Database Key:
463: . -pc_mg_levels levels - set the number of levels to use

465:   Level: intermediate

467:   Notes:
468:   If the number of levels is one then the multigrid uses the `-mg_levels` prefix
469:   for setting the level options rather than the `-mg_coarse` or `-mg_fine` prefix.

471:   You can free the information in comms after this routine is called.

473:   The array of MPI communicators must contain `MPI_COMM_NULL` for those ranks that at each level
474:   are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on
475:   the two levels, rank 0 in the original communicator will pass in an array of 2 communicators
476:   of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators
477:   the first of size 2 and the second of value `MPI_COMM_NULL` since the rank 1 does not participate
478:   in the coarse grid solve.

480:   Since each coarser level may have a new `MPI_Comm` with fewer ranks than the previous, one
481:   must take special care in providing the restriction and interpolation operation. We recommend
482:   providing these as two step operations; first perform a standard restriction or interpolation on
483:   the full number of ranks for that level and then use an MPI call to copy the resulting vector
484:   array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, note in both
485:   cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and
486:   receives or `MPI_AlltoAllv()` could be used to do the reshuffling of the vector entries.

488:   Fortran Notes:
489:   Use comms = `PETSC_NULL_MPI_COMM` as the equivalent of `NULL` in the C interface. Note `PETSC_NULL_MPI_COMM`
490:   is not `MPI_COMM_NULL`. It is more like `PETSC_NULL_INTEGER`, `PETSC_NULL_REAL` etc.

492: .seealso: [](ch_ksp), `PCMGSetType()`, `PCMGGetLevels()`
493: @*/
494: PetscErrorCode PCMGSetLevels(PC pc, PetscInt levels, MPI_Comm *comms)
495: {
496:   PetscFunctionBegin;
498:   if (comms) PetscAssertPointer(comms, 3);
499:   PetscTryMethod(pc, "PCMGSetLevels_C", (PC, PetscInt, MPI_Comm *), (pc, levels, comms));
500:   PetscFunctionReturn(PETSC_SUCCESS);
501: }

503: PetscErrorCode PCDestroy_MG(PC pc)
504: {
505:   PC_MG         *mg       = (PC_MG *)pc->data;
506:   PC_MG_Levels **mglevels = mg->levels;
507:   PetscInt       i, n;

509:   PetscFunctionBegin;
510:   PetscCall(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) PetscCall(KSPDestroy(&mglevels[i]->smoothd));
515:       PetscCall(KSPDestroy(&mglevels[i]->smoothu));
516:       PetscCall(KSPDestroy(&mglevels[i]->cr));
517:       PetscCall(PetscFree(mglevels[i]));
518:     }
519:     PetscCall(PetscFree(mg->levels));
520:   }
521:   PetscCall(PetscFree(pc->data));
522:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL));
523:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL));
524:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", NULL));
525:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetReusePreconditioner_C", NULL));
526:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", NULL));
527:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", NULL));
528:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL));
529:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL));
530:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", NULL));
531:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", NULL));
532:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", NULL));
533:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", NULL));
534:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", NULL));
535:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", NULL));
536:   PetscFunctionReturn(PETSC_SUCCESS);
537: }

539: /*
540:    PCApply_MG - Runs either an additive, multiplicative, Kaskadic
541:              or full cycle of multigrid.

543:   Note:
544:   A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle().
545: */
546: static PetscErrorCode PCApply_MG_Internal(PC pc, Vec b, Vec x, Mat B, Mat X, PetscBool transpose)
547: {
548:   PC_MG         *mg       = (PC_MG *)pc->data;
549:   PC_MG_Levels **mglevels = mg->levels;
550:   PC             tpc;
551:   PetscInt       levels = mglevels[0]->levels, i;
552:   PetscBool      changeu, changed, matapp;

554:   PetscFunctionBegin;
555:   matapp = (PetscBool)(B && X);
556:   if (mg->stageApply) PetscCall(PetscLogStagePush(mg->stageApply));
557:   /* When the DM is supplying the matrix then it will not exist until here */
558:   for (i = 0; i < levels; i++) {
559:     if (!mglevels[i]->A) {
560:       PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL));
561:       PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A));
562:     }
563:   }

565:   PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc));
566:   PetscCall(PCPreSolveChangeRHS(tpc, &changed));
567:   PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc));
568:   PetscCall(PCPreSolveChangeRHS(tpc, &changeu));
569:   if (!changeu && !changed) {
570:     if (matapp) {
571:       PetscCall(MatDestroy(&mglevels[levels - 1]->B));
572:       mglevels[levels - 1]->B = B;
573:     } else {
574:       PetscCall(VecDestroy(&mglevels[levels - 1]->b));
575:       mglevels[levels - 1]->b = b;
576:     }
577:   } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */
578:     if (matapp) {
579:       if (mglevels[levels - 1]->B) {
580:         PetscInt  N1, N2;
581:         PetscBool flg;

583:         PetscCall(MatGetSize(mglevels[levels - 1]->B, NULL, &N1));
584:         PetscCall(MatGetSize(B, NULL, &N2));
585:         PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels - 1]->B, ((PetscObject)B)->type_name, &flg));
586:         if (N1 != N2 || !flg) PetscCall(MatDestroy(&mglevels[levels - 1]->B));
587:       }
588:       if (!mglevels[levels - 1]->B) PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &mglevels[levels - 1]->B));
589:       else PetscCall(MatCopy(B, mglevels[levels - 1]->B, SAME_NONZERO_PATTERN));
590:     } else {
591:       if (!mglevels[levels - 1]->b) {
592:         Vec *vec;

594:         PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL));
595:         mglevels[levels - 1]->b = *vec;
596:         PetscCall(PetscFree(vec));
597:       }
598:       PetscCall(VecCopy(b, mglevels[levels - 1]->b));
599:     }
600:   }
601:   if (matapp) {
602:     mglevels[levels - 1]->X = X;
603:   } else {
604:     mglevels[levels - 1]->x = x;
605:   }

607:   /* If coarser Xs are present, it means we have already block applied the PC at least once
608:      Reset operators if sizes/type do no match */
609:   if (matapp && levels > 1 && mglevels[levels - 2]->X) {
610:     PetscInt  Xc, Bc;
611:     PetscBool flg;

613:     PetscCall(MatGetSize(mglevels[levels - 2]->X, NULL, &Xc));
614:     PetscCall(MatGetSize(mglevels[levels - 1]->B, NULL, &Bc));
615:     PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels - 2]->X, ((PetscObject)mglevels[levels - 1]->X)->type_name, &flg));
616:     if (Xc != Bc || !flg) {
617:       PetscCall(MatDestroy(&mglevels[levels - 1]->R));
618:       for (i = 0; i < levels - 1; i++) {
619:         PetscCall(MatDestroy(&mglevels[i]->R));
620:         PetscCall(MatDestroy(&mglevels[i]->B));
621:         PetscCall(MatDestroy(&mglevels[i]->X));
622:       }
623:     }
624:   }

626:   if (mg->am == PC_MG_MULTIPLICATIVE) {
627:     if (matapp) PetscCall(MatZeroEntries(X));
628:     else PetscCall(VecZeroEntries(x));
629:     for (i = 0; i < mg->cyclesperpcapply; i++) PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, transpose, matapp, NULL));
630:   } else if (mg->am == PC_MG_ADDITIVE) {
631:     PetscCall(PCMGACycle_Private(pc, mglevels, transpose, matapp));
632:   } else if (mg->am == PC_MG_KASKADE) {
633:     PetscCall(PCMGKCycle_Private(pc, mglevels, transpose, matapp));
634:   } else {
635:     PetscCall(PCMGFCycle_Private(pc, mglevels, transpose, matapp));
636:   }
637:   if (mg->stageApply) PetscCall(PetscLogStagePop());
638:   if (!changeu && !changed) {
639:     if (matapp) {
640:       mglevels[levels - 1]->B = NULL;
641:     } else {
642:       mglevels[levels - 1]->b = NULL;
643:     }
644:   }
645:   PetscFunctionReturn(PETSC_SUCCESS);
646: }

648: static PetscErrorCode PCApply_MG(PC pc, Vec b, Vec x)
649: {
650:   PetscFunctionBegin;
651:   PetscCall(PCApply_MG_Internal(pc, b, x, NULL, NULL, PETSC_FALSE));
652:   PetscFunctionReturn(PETSC_SUCCESS);
653: }

655: static PetscErrorCode PCApplyTranspose_MG(PC pc, Vec b, Vec x)
656: {
657:   PetscFunctionBegin;
658:   PetscCall(PCApply_MG_Internal(pc, b, x, NULL, NULL, PETSC_TRUE));
659:   PetscFunctionReturn(PETSC_SUCCESS);
660: }

662: static PetscErrorCode PCMatApply_MG(PC pc, Mat b, Mat x)
663: {
664:   PetscFunctionBegin;
665:   PetscCall(PCApply_MG_Internal(pc, NULL, NULL, b, x, PETSC_FALSE));
666:   PetscFunctionReturn(PETSC_SUCCESS);
667: }

669: static PetscErrorCode PCMatApplyTranspose_MG(PC pc, Mat b, Mat x)
670: {
671:   PetscFunctionBegin;
672:   PetscCall(PCApply_MG_Internal(pc, NULL, NULL, b, x, PETSC_TRUE));
673:   PetscFunctionReturn(PETSC_SUCCESS);
674: }

676: PetscErrorCode PCSetFromOptions_MG(PC pc, PetscOptionItems PetscOptionsObject)
677: {
678:   PetscInt            levels, cycles;
679:   PetscBool           flg, flg2;
680:   PC_MG              *mg = (PC_MG *)pc->data;
681:   PC_MG_Levels      **mglevels;
682:   PCMGType            mgtype;
683:   PCMGCycleType       mgctype;
684:   PCMGGalerkinType    gtype;
685:   PCMGCoarseSpaceType coarseSpaceType;

687:   PetscFunctionBegin;
688:   levels = PetscMax(mg->nlevels, 1);
689:   PetscOptionsHeadBegin(PetscOptionsObject, "Multigrid options");
690:   PetscCall(PetscOptionsInt("-pc_mg_levels", "Number of Levels", "PCMGSetLevels", levels, &levels, &flg));
691:   if (!flg && !mg->levels && pc->dm) {
692:     PetscCall(DMGetRefineLevel(pc->dm, &levels));
693:     levels++;
694:     mg->usedmfornumberoflevels = PETSC_TRUE;
695:   }
696:   PetscCall(PCMGSetLevels(pc, levels, NULL));
697:   mglevels = mg->levels;

699:   mgctype = (PCMGCycleType)mglevels[0]->cycles;
700:   PetscCall(PetscOptionsEnum("-pc_mg_cycle_type", "V cycle or for W-cycle", "PCMGSetCycleType", PCMGCycleTypes, (PetscEnum)mgctype, (PetscEnum *)&mgctype, &flg));
701:   if (flg) PetscCall(PCMGSetCycleType(pc, mgctype));
702:   coarseSpaceType = mg->coarseSpaceType;
703:   PetscCall(PetscOptionsEnum("-pc_mg_adapt_interp_coarse_space", "Type of adaptive coarse space: none, polynomial, harmonic, eigenvector, generalized_eigenvector, gdsw", "PCMGSetAdaptCoarseSpaceType", PCMGCoarseSpaceTypes, (PetscEnum)coarseSpaceType, (PetscEnum *)&coarseSpaceType, &flg));
704:   if (flg) PetscCall(PCMGSetAdaptCoarseSpaceType(pc, coarseSpaceType));
705:   PetscCall(PetscOptionsInt("-pc_mg_adapt_interp_n", "Size of the coarse space for adaptive interpolation", "PCMGSetCoarseSpace", mg->Nc, &mg->Nc, &flg));
706:   PetscCall(PetscOptionsBool("-pc_mg_mesp_monitor", "Monitor the multilevel eigensolver", "PCMGSetAdaptInterpolation", PETSC_FALSE, &mg->mespMonitor, &flg));
707:   flg2 = PETSC_FALSE;
708:   PetscCall(PetscOptionsBool("-pc_mg_adapt_cr", "Monitor coarse space quality using Compatible Relaxation (CR)", "PCMGSetAdaptCR", PETSC_FALSE, &flg2, &flg));
709:   if (flg) PetscCall(PCMGSetAdaptCR(pc, flg2));
710:   flg = PETSC_FALSE;
711:   PetscCall(PetscOptionsBool("-pc_mg_distinct_smoothup", "Create separate smoothup KSP and append the prefix _up", "PCMGSetDistinctSmoothUp", PETSC_FALSE, &flg, NULL));
712:   if (flg) PetscCall(PCMGSetDistinctSmoothUp(pc));
713:   PetscCall(PetscOptionsEnum("-pc_mg_galerkin", "Use Galerkin process to compute coarser operators", "PCMGSetGalerkin", PCMGGalerkinTypes, (PetscEnum)mg->galerkin, (PetscEnum *)&gtype, &flg));
714:   if (flg) PetscCall(PCMGSetGalerkin(pc, gtype));
715:   mgtype = mg->am;
716:   PetscCall(PetscOptionsEnum("-pc_mg_type", "Multigrid type", "PCMGSetType", PCMGTypes, (PetscEnum)mgtype, (PetscEnum *)&mgtype, &flg));
717:   if (flg) PetscCall(PCMGSetType(pc, mgtype));
718:   if (mg->am == PC_MG_MULTIPLICATIVE) {
719:     PetscCall(PetscOptionsInt("-pc_mg_multiplicative_cycles", "Number of cycles for each preconditioner step", "PCMGMultiplicativeSetCycles", mg->cyclesperpcapply, &cycles, &flg));
720:     if (flg) PetscCall(PCMGMultiplicativeSetCycles(pc, cycles));
721:   }
722:   flg = PETSC_FALSE;
723:   PetscCall(PetscOptionsBool("-pc_mg_log", "Log times for each multigrid level", "None", flg, &flg, NULL));
724:   if (flg) {
725:     PetscInt i;
726:     char     eventname[128];

728:     levels = mglevels[0]->levels;
729:     for (i = 0; i < levels; i++) {
730:       PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGSetup Level %" PetscInt_FMT, i));
731:       PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventsmoothsetup));
732:       PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGSmooth Level %" PetscInt_FMT, i));
733:       PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventsmoothsolve));
734:       if (i) {
735:         PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGResid Level %" PetscInt_FMT, i));
736:         PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventresidual));
737:         PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGInterp Level %" PetscInt_FMT, i));
738:         PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventinterprestrict));
739:       }
740:     }

742:     if (PetscDefined(USE_LOG)) {
743:       const char sname[] = "MG Apply";

745:       PetscCall(PetscLogStageGetId(sname, &mg->stageApply));
746:       if (mg->stageApply < 0) PetscCall(PetscLogStageRegister(sname, &mg->stageApply));
747:     }
748:   }
749:   PetscOptionsHeadEnd();
750:   /* Check option consistency */
751:   PetscCall(PCMGGetGalerkin(pc, &gtype));
752:   PetscCall(PCMGGetAdaptInterpolation(pc, &flg));
753:   PetscCheck(!flg || !(gtype >= PC_MG_GALERKIN_NONE), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator");
754:   PetscFunctionReturn(PETSC_SUCCESS);
755: }

757: const char *const PCMGTypes[]            = {"MULTIPLICATIVE", "ADDITIVE", "FULL", "KASKADE", "PCMGType", "PC_MG", NULL};
758: const char *const PCMGCycleTypes[]       = {"invalid", "v", "w", "PCMGCycleType", "PC_MG_CYCLE", NULL};
759: const char *const PCMGGalerkinTypes[]    = {"both", "pmat", "mat", "none", "external", "PCMGGalerkinType", "PC_MG_GALERKIN", NULL};
760: const char *const PCMGCoarseSpaceTypes[] = {"none", "polynomial", "harmonic", "eigenvector", "generalized_eigenvector", "gdsw", "PCMGCoarseSpaceType", "PCMG_ADAPT_NONE", NULL};

762: #include <petscdraw.h>
763: PetscErrorCode PCView_MG(PC pc, PetscViewer viewer)
764: {
765:   PC_MG         *mg       = (PC_MG *)pc->data;
766:   PC_MG_Levels **mglevels = mg->levels;
767:   PetscInt       levels   = mglevels ? mglevels[0]->levels : 0, i;
768:   PetscBool      isascii, isbinary, isdraw;

770:   PetscFunctionBegin;
771:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
772:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
773:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
774:   if (isascii) {
775:     const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown";

777:     if (levels == 1) PetscCall(PetscViewerASCIIPrintf(viewer, "  WARNING: Multigrid is being run with only a single level!\n"));
778:     PetscCall(PetscViewerASCIIPrintf(viewer, "  type is %s, levels=%" PetscInt_FMT " cycles=%s\n", PCMGTypes[mg->am], levels, cyclename));
779:     if (mg->am == PC_MG_MULTIPLICATIVE) PetscCall(PetscViewerASCIIPrintf(viewer, "    Cycles per PCApply=%" PetscInt_FMT "\n", mg->cyclesperpcapply));
780:     if (mg->galerkin == PC_MG_GALERKIN_BOTH) {
781:       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using Galerkin computed coarse grid matrices\n"));
782:     } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) {
783:       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using Galerkin computed coarse grid matrices for pmat\n"));
784:     } else if (mg->galerkin == PC_MG_GALERKIN_MAT) {
785:       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using Galerkin computed coarse grid matrices for mat\n"));
786:     } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) {
787:       PetscCall(PetscViewerASCIIPrintf(viewer, "    Using externally compute Galerkin coarse grid matrices\n"));
788:     } else {
789:       PetscCall(PetscViewerASCIIPrintf(viewer, "    Not using Galerkin computed coarse grid matrices\n"));
790:     }
791:     if (mg->view) PetscCall((*mg->view)(pc, viewer));
792:     for (i = 0; i < levels; i++) {
793:       if (i) {
794:         PetscCall(PetscViewerASCIIPrintf(viewer, "Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n", i));
795:       } else {
796:         PetscCall(PetscViewerASCIIPrintf(viewer, "Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n", i));
797:       }
798:       PetscCall(PetscViewerASCIIPushTab(viewer));
799:       PetscCall(KSPView(mglevels[i]->smoothd, viewer));
800:       PetscCall(PetscViewerASCIIPopTab(viewer));
801:       if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) {
802:         PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) same as down solver (pre-smoother)\n"));
803:       } else if (i) {
804:         PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n", i));
805:         PetscCall(PetscViewerASCIIPushTab(viewer));
806:         PetscCall(KSPView(mglevels[i]->smoothu, viewer));
807:         PetscCall(PetscViewerASCIIPopTab(viewer));
808:       }
809:       if (i && mglevels[i]->cr) {
810:         PetscCall(PetscViewerASCIIPrintf(viewer, "CR solver on level %" PetscInt_FMT " -------------------------------\n", i));
811:         PetscCall(PetscViewerASCIIPushTab(viewer));
812:         PetscCall(KSPView(mglevels[i]->cr, viewer));
813:         PetscCall(PetscViewerASCIIPopTab(viewer));
814:       }
815:     }
816:   } else if (isbinary) {
817:     for (i = levels - 1; i >= 0; i--) {
818:       PetscCall(KSPView(mglevels[i]->smoothd, viewer));
819:       if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPView(mglevels[i]->smoothu, viewer));
820:     }
821:   } else if (isdraw) {
822:     PetscDraw draw;
823:     PetscReal x, w, y, bottom, th;
824:     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
825:     PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y));
826:     PetscCall(PetscDrawStringGetSize(draw, NULL, &th));
827:     bottom = y - th;
828:     for (i = levels - 1; i >= 0; i--) {
829:       if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) {
830:         PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom));
831:         PetscCall(KSPView(mglevels[i]->smoothd, viewer));
832:         PetscCall(PetscDrawPopCurrentPoint(draw));
833:       } else {
834:         w = 0.5 * PetscMin(1.0 - x, x);
835:         PetscCall(PetscDrawPushCurrentPoint(draw, x + w, bottom));
836:         PetscCall(KSPView(mglevels[i]->smoothd, viewer));
837:         PetscCall(PetscDrawPopCurrentPoint(draw));
838:         PetscCall(PetscDrawPushCurrentPoint(draw, x - w, bottom));
839:         PetscCall(KSPView(mglevels[i]->smoothu, viewer));
840:         PetscCall(PetscDrawPopCurrentPoint(draw));
841:       }
842:       PetscCall(PetscDrawGetBoundingBox(draw, NULL, &bottom, NULL, NULL));
843:       bottom -= th;
844:     }
845:   }
846:   PetscFunctionReturn(PETSC_SUCCESS);
847: }

849: #include <petsc/private/kspimpl.h>

851: /*
852:     Calls setup for the KSP on each level
853: */
854: PetscErrorCode PCSetUp_MG(PC pc)
855: {
856:   PC_MG         *mg       = (PC_MG *)pc->data;
857:   PC_MG_Levels **mglevels = mg->levels;
858:   PetscInt       n;
859:   PC             cpc;
860:   PetscBool      dump = PETSC_FALSE, opsset, use_amat, missinginterpolate = PETSC_FALSE;
861:   Mat            dA, dB;
862:   Vec            tvec;
863:   DM            *dms;
864:   PetscViewer    viewer = NULL;
865:   PetscBool      dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE;
866:   PetscBool      adaptInterpolation = mg->adaptInterpolation;

868:   PetscFunctionBegin;
869:   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels with PCMGSetLevels() before setting up");
870:   n = mglevels[0]->levels;
871:   /* FIX: Move this to PCSetFromOptions_MG? */
872:   if (mg->usedmfornumberoflevels) {
873:     PetscInt levels;

875:     PetscCall(DMGetRefineLevel(pc->dm, &levels));
876:     levels++;
877:     if (levels > n) { /* the problem is now being solved on a finer grid */
878:       PetscCall(PCMGSetLevels(pc, levels, NULL));
879:       n = levels;
880:       PetscCall(PCSetFromOptions(pc)); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */
881:       mglevels = mg->levels;
882:     }
883:   }

885:   /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */
886:   /* so use those from global PC */
887:   /* Is this what we always want? What if user wants to keep old one? */
888:   PetscCall(KSPGetOperatorsSet(mglevels[n - 1]->smoothd, NULL, &opsset));
889:   if (opsset) {
890:     Mat mmat;
891:     PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, NULL, &mmat));
892:     if (mmat == pc->pmat) opsset = PETSC_FALSE;
893:   }
894:   /* fine grid smoother inherits the reuse-pc flag */
895:   PetscCall(KSPGetPC(mglevels[n - 1]->smoothd, &cpc));
896:   cpc->reusepreconditioner = pc->reusepreconditioner;
897:   PetscCall(KSPGetPC(mglevels[n - 1]->smoothu, &cpc));
898:   cpc->reusepreconditioner = pc->reusepreconditioner;

900:   /* Create CR solvers */
901:   PetscCall(PCMGGetAdaptCR(pc, &doCR));
902:   if (doCR) {
903:     const char *prefix;

905:     PetscCall(PCGetOptionsPrefix(pc, &prefix));
906:     for (PetscInt i = 1; i < n; ++i) {
907:       PC   ipc, cr;
908:       char crprefix[128];

910:       PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &mglevels[i]->cr));
911:       PetscCall(KSPSetNestLevel(mglevels[i]->cr, pc->kspnestlevel));
912:       PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE));
913:       PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->cr, (PetscObject)pc, n - i));
914:       PetscCall(KSPSetOptionsPrefix(mglevels[i]->cr, prefix));
915:       PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level));
916:       PetscCall(KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV));
917:       PetscCall(KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL));
918:       PetscCall(KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED));
919:       PetscCall(KSPGetPC(mglevels[i]->cr, &ipc));

921:       PetscCall(PCSetType(ipc, PCCOMPOSITE));
922:       PetscCall(PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE));
923:       PetscCall(PCCompositeAddPCType(ipc, PCSOR));
924:       PetscCall(CreateCR_Private(pc, i, &cr));
925:       PetscCall(PCCompositeAddPC(ipc, cr));
926:       PetscCall(PCDestroy(&cr));

928:       PetscCall(KSPSetTolerances(mglevels[i]->cr, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, mg->default_smoothd));
929:       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE));
930:       PetscCall(PetscSNPrintf(crprefix, 128, "mg_levels_%" PetscInt_FMT "_cr_", i));
931:       PetscCall(KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix));
932:     }
933:   }

935:   if (!opsset) {
936:     PetscCall(PCGetUseAmat(pc, &use_amat));
937:     if (use_amat) {
938:       PetscCall(PetscInfo(pc, "Using outer operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n"));
939:       PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->mat, pc->pmat));
940:     } else {
941:       PetscCall(PetscInfo(pc, "Using matrix (pmat) operators to define finest grid operator \n  because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n"));
942:       PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->pmat, pc->pmat));
943:     }
944:   }

946:   for (PetscInt i = n - 1; i > 0; i--) {
947:     if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) {
948:       missinginterpolate = PETSC_TRUE;
949:       break;
950:     }
951:   }

953:   PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, &dA, &dB));
954:   if (dA == dB) dAeqdB = PETSC_TRUE;
955:   if (mg->galerkin == PC_MG_GALERKIN_NONE || ((mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_MAT) && !dAeqdB)) needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */

957:   if (pc->dm && !pc->setupcalled) {
958:     /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */
959:     PetscCall(KSPSetDM(mglevels[n - 1]->smoothd, pc->dm));
960:     PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothd, KSP_DMACTIVE_ALL, PETSC_FALSE));
961:     if (mglevels[n - 1]->smoothd != mglevels[n - 1]->smoothu) {
962:       PetscCall(KSPSetDM(mglevels[n - 1]->smoothu, pc->dm));
963:       PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothu, KSP_DMACTIVE_ALL, PETSC_FALSE));
964:     }
965:     if (mglevels[n - 1]->cr) {
966:       PetscCall(KSPSetDM(mglevels[n - 1]->cr, pc->dm));
967:       PetscCall(KSPSetDMActive(mglevels[n - 1]->cr, KSP_DMACTIVE_ALL, PETSC_FALSE));
968:     }
969:   }

971:   /*
972:    Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS)
973:    Skipping for externally managed hierarchy (such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs?
974:   */
975:   if (missinginterpolate && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) {
976:     /* first see if we can compute a coarse space */
977:     if (mg->coarseSpaceType == PCMG_ADAPT_GDSW) {
978:       for (PetscInt i = n - 2; i > -1; i--) {
979:         if (!mglevels[i + 1]->restrct && !mglevels[i + 1]->interpolate) {
980:           PetscCall(PCMGComputeCoarseSpace_Internal(pc, i + 1, mg->coarseSpaceType, mg->Nc, NULL, &mglevels[i + 1]->coarseSpace));
981:           PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->coarseSpace));
982:         }
983:       }
984:     } else { /* construct the interpolation from the DMs */
985:       Mat p;
986:       Vec rscale;

988:       PetscCheck(n == 1 || pc->dm, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "PC lacks a DM so cannot automatically construct a multigrid hierarchy. Number of levels requested %" PetscInt_FMT, n);
989:       PetscCall(PetscMalloc1(n, &dms));
990:       dms[n - 1] = pc->dm;
991:       /* Separately create them so we do not get DMKSP interference between levels */
992:       for (PetscInt i = n - 2; i > -1; i--) PetscCall(DMCoarsen(dms[i + 1], MPI_COMM_NULL, &dms[i]));
993:       for (PetscInt i = n - 2; i > -1; i--) {
994:         PetscBool dmhasrestrict, dmhasinject;

996:         PetscCall(KSPSetDM(mglevels[i]->smoothd, dms[i]));
997:         if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothd, KSP_DMACTIVE_ALL, PETSC_FALSE));
998:         PetscCall(KSPSetDMActive(mglevels[i]->smoothd, KSP_DMACTIVE_RHS, PETSC_FALSE));
999:         if (mglevels[i]->smoothd != mglevels[i]->smoothu) {
1000:           PetscCall(KSPSetDM(mglevels[i]->smoothu, dms[i]));
1001:           if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothu, KSP_DMACTIVE_ALL, PETSC_FALSE));
1002:           PetscCall(KSPSetDMActive(mglevels[i]->smoothu, KSP_DMACTIVE_RHS, PETSC_FALSE));
1003:         }
1004:         if (mglevels[i]->cr) {
1005:           PetscCall(KSPSetDM(mglevels[i]->cr, dms[i]));
1006:           if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->cr, KSP_DMACTIVE_ALL, PETSC_FALSE));
1007:           PetscCall(KSPSetDMActive(mglevels[i]->cr, KSP_DMACTIVE_RHS, PETSC_FALSE));
1008:         }
1009:         if (!mglevels[i + 1]->interpolate) {
1010:           PetscCall(DMCreateInterpolation(dms[i], dms[i + 1], &p, &rscale));
1011:           PetscCall(PCMGSetInterpolation(pc, i + 1, p));
1012:           if (rscale) PetscCall(PCMGSetRScale(pc, i + 1, rscale));
1013:           PetscCall(VecDestroy(&rscale));
1014:           PetscCall(MatDestroy(&p));
1015:         }
1016:         PetscCall(DMHasCreateRestriction(dms[i], &dmhasrestrict));
1017:         if (dmhasrestrict && !mglevels[i + 1]->restrct) {
1018:           PetscCall(DMCreateRestriction(dms[i], dms[i + 1], &p));
1019:           PetscCall(PCMGSetRestriction(pc, i + 1, p));
1020:           PetscCall(MatDestroy(&p));
1021:         }
1022:         PetscCall(DMHasCreateInjection(dms[i], &dmhasinject));
1023:         if (dmhasinject && !mglevels[i + 1]->inject) {
1024:           PetscCall(DMCreateInjection(dms[i], dms[i + 1], &p));
1025:           PetscCall(PCMGSetInjection(pc, i + 1, p));
1026:           PetscCall(MatDestroy(&p));
1027:         }
1028:       }

1030:       for (PetscInt i = n - 2; i > -1; i--) PetscCall(DMDestroy(&dms[i]));
1031:       PetscCall(PetscFree(dms));
1032:     }
1033:   }

1035:   if (mg->galerkin < PC_MG_GALERKIN_NONE) {
1036:     Mat       A, B;
1037:     PetscBool doA = PETSC_FALSE, doB = PETSC_FALSE;
1038:     MatReuse  reuse = MAT_INITIAL_MATRIX;

1040:     if (mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_BOTH) doB = PETSC_TRUE;
1041:     if (mg->galerkin == PC_MG_GALERKIN_MAT || (mg->galerkin == PC_MG_GALERKIN_BOTH && dA != dB)) doA = PETSC_TRUE;
1042:     if (pc->setupcalled) reuse = MAT_REUSE_MATRIX;
1043:     for (PetscInt i = n - 2; i > -1; i--) {
1044:       PetscCheck(mglevels[i + 1]->restrct || mglevels[i + 1]->interpolate, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must provide interpolation or restriction for each MG level except level 0");
1045:       if (!mglevels[i + 1]->interpolate) PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->restrct));
1046:       if (!mglevels[i + 1]->restrct) PetscCall(PCMGSetRestriction(pc, i + 1, mglevels[i + 1]->interpolate));
1047:       if (reuse == MAT_REUSE_MATRIX) PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, &B));
1048:       if (doA) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dA, mglevels[i + 1]->interpolate, reuse, 1.0, &A));
1049:       if (doB) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dB, mglevels[i + 1]->interpolate, reuse, 1.0, &B));
1050:       /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */
1051:       if (!doA && dAeqdB) {
1052:         if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)B));
1053:         A = B;
1054:       } else if (!doA && reuse == MAT_INITIAL_MATRIX) {
1055:         PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, NULL));
1056:         PetscCall(PetscObjectReference((PetscObject)A));
1057:       }
1058:       if (!doB && dAeqdB) {
1059:         if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)A));
1060:         B = A;
1061:       } else if (!doB && reuse == MAT_INITIAL_MATRIX) {
1062:         PetscCall(KSPGetOperators(mglevels[i]->smoothd, NULL, &B));
1063:         PetscCall(PetscObjectReference((PetscObject)B));
1064:       }
1065:       if (reuse == MAT_INITIAL_MATRIX) {
1066:         PetscCall(KSPSetOperators(mglevels[i]->smoothd, A, B));
1067:         PetscCall(PetscObjectDereference((PetscObject)A));
1068:         PetscCall(PetscObjectDereference((PetscObject)B));
1069:       }
1070:       dA = A;
1071:       dB = B;
1072:     }
1073:   }

1075:   /* Adapt interpolation matrices */
1076:   if (adaptInterpolation) {
1077:     for (PetscInt i = 0; i < n; ++i) {
1078:       if (!mglevels[i]->coarseSpace) PetscCall(PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i - 1]->coarseSpace, &mglevels[i]->coarseSpace));
1079:       if (i) PetscCall(PCMGAdaptInterpolator_Internal(pc, i, mglevels[i - 1]->smoothu, mglevels[i]->smoothu, mglevels[i - 1]->coarseSpace, mglevels[i]->coarseSpace));
1080:     }
1081:     for (PetscInt i = n - 2; i > -1; --i) PetscCall(PCMGRecomputeLevelOperators_Internal(pc, i));
1082:   }

1084:   if (needRestricts && pc->dm) {
1085:     for (PetscInt i = n - 2; i >= 0; i--) {
1086:       DM  dmfine, dmcoarse;
1087:       Mat Restrict, Inject;
1088:       Vec rscale;

1090:       PetscCall(KSPGetDM(mglevels[i + 1]->smoothd, &dmfine));
1091:       PetscCall(KSPGetDM(mglevels[i]->smoothd, &dmcoarse));
1092:       PetscCall(PCMGGetRestriction(pc, i + 1, &Restrict));
1093:       PetscCall(PCMGGetRScale(pc, i + 1, &rscale));
1094:       PetscCall(PCMGGetInjection(pc, i + 1, &Inject));
1095:       PetscCall(DMRestrict(dmfine, Restrict, rscale, Inject, dmcoarse));
1096:     }
1097:   }

1099:   if (!pc->setupcalled) {
1100:     for (PetscInt i = 0; i < n; i++) PetscCall(KSPSetFromOptions(mglevels[i]->smoothd));
1101:     for (PetscInt i = 1; i < n; i++) {
1102:       if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) PetscCall(KSPSetFromOptions(mglevels[i]->smoothu));
1103:       if (mglevels[i]->cr) PetscCall(KSPSetFromOptions(mglevels[i]->cr));
1104:     }
1105:     /* insure that if either interpolation or restriction is set the other one is set */
1106:     for (PetscInt i = 1; i < n; i++) {
1107:       PetscCall(PCMGGetInterpolation(pc, i, NULL));
1108:       PetscCall(PCMGGetRestriction(pc, i, NULL));
1109:     }
1110:     for (PetscInt i = 0; i < n - 1; i++) {
1111:       if (!mglevels[i]->b) {
1112:         Vec *vec;
1113:         PetscCall(KSPCreateVecs(mglevels[i]->smoothd, 1, &vec, 0, NULL));
1114:         PetscCall(PCMGSetRhs(pc, i, *vec));
1115:         PetscCall(VecDestroy(vec));
1116:         PetscCall(PetscFree(vec));
1117:       }
1118:       if (!mglevels[i]->r && i) {
1119:         PetscCall(VecDuplicate(mglevels[i]->b, &tvec));
1120:         PetscCall(PCMGSetR(pc, i, tvec));
1121:         PetscCall(VecDestroy(&tvec));
1122:       }
1123:       if (!mglevels[i]->x) {
1124:         PetscCall(VecDuplicate(mglevels[i]->b, &tvec));
1125:         PetscCall(PCMGSetX(pc, i, tvec));
1126:         PetscCall(VecDestroy(&tvec));
1127:       }
1128:       if (doCR) {
1129:         PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crx));
1130:         PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crb));
1131:       }
1132:     }
1133:     if (n != 1 && !mglevels[n - 1]->r) {
1134:       /* PCMGSetR() on the finest level if user did not supply it */
1135:       Vec *vec;

1137:       PetscCall(KSPCreateVecs(mglevels[n - 1]->smoothd, 1, &vec, 0, NULL));
1138:       PetscCall(PCMGSetR(pc, n - 1, *vec));
1139:       PetscCall(VecDestroy(vec));
1140:       PetscCall(PetscFree(vec));
1141:     }
1142:     if (doCR) {
1143:       PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crx));
1144:       PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crb));
1145:     }
1146:   }

1148:   if (pc->dm) {
1149:     /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */
1150:     for (PetscInt i = 0; i < n - 1; i++) {
1151:       if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;
1152:     }
1153:   }
1154:   // We got here (PCSetUp_MG) because the matrix has changed, which means the smoother needs to be set up again (e.g.,
1155:   // new diagonal for Jacobi). Setting it here allows it to be logged under PCSetUp rather than deep inside a PCApply.
1156:   if (mglevels[n - 1]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[n - 1]->smoothd->setupstage = KSP_SETUP_NEWMATRIX;

1158:   for (PetscInt i = 1; i < n; i++) {
1159:     if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1) {
1160:       /* if doing only down then initial guess is zero */
1161:       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE));
1162:     }
1163:     if (mglevels[i]->cr) PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE));
1164:     if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
1165:     PetscCall(KSPSetUp(mglevels[i]->smoothd));
1166:     if (mglevels[i]->smoothd->reason) pc->failedreason = PC_SUBPC_ERROR;
1167:     if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
1168:     if (!mglevels[i]->residual) {
1169:       Mat mat;

1171:       PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL));
1172:       PetscCall(PCMGSetResidual(pc, i, PCMGResidualDefault, mat));
1173:     }
1174:     if (!mglevels[i]->residualtranspose) {
1175:       Mat mat;

1177:       PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL));
1178:       PetscCall(PCMGSetResidualTranspose(pc, i, PCMGResidualTransposeDefault, mat));
1179:     }
1180:   }
1181:   for (PetscInt i = 1; i < n; i++) {
1182:     if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) {
1183:       Mat downmat, downpmat;

1185:       /* check if operators have been set for up, if not use down operators to set them */
1186:       PetscCall(KSPGetOperatorsSet(mglevels[i]->smoothu, &opsset, NULL));
1187:       if (!opsset) {
1188:         PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat));
1189:         PetscCall(KSPSetOperators(mglevels[i]->smoothu, downmat, downpmat));
1190:       }

1192:       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothu, PETSC_TRUE));
1193:       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
1194:       PetscCall(KSPSetUp(mglevels[i]->smoothu));
1195:       if (mglevels[i]->smoothu->reason) pc->failedreason = PC_SUBPC_ERROR;
1196:       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
1197:     }
1198:     if (mglevels[i]->cr) {
1199:       Mat downmat, downpmat;

1201:       /* check if operators have been set for up, if not use down operators to set them */
1202:       PetscCall(KSPGetOperatorsSet(mglevels[i]->cr, &opsset, NULL));
1203:       if (!opsset) {
1204:         PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat));
1205:         PetscCall(KSPSetOperators(mglevels[i]->cr, downmat, downpmat));
1206:       }

1208:       PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE));
1209:       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
1210:       PetscCall(KSPSetUp(mglevels[i]->cr));
1211:       if (mglevels[i]->cr->reason) pc->failedreason = PC_SUBPC_ERROR;
1212:       if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0));
1213:     }
1214:   }

1216:   if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0));
1217:   PetscCall(KSPSetUp(mglevels[0]->smoothd));
1218:   if (mglevels[0]->smoothd->reason) pc->failedreason = PC_SUBPC_ERROR;
1219:   if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0));

1221:   /*
1222:      Dump the interpolation/restriction matrices plus the
1223:    Jacobian/stiffness on each level. This allows MATLAB users to
1224:    easily check if the Galerkin condition A_c = R A_f R^T is satisfied.

1226:    Only support one or the other at the same time.
1227:   */
1228: #if defined(PETSC_USE_SOCKET_VIEWER)
1229:   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_matlab", &dump, NULL));
1230:   if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc));
1231:   dump = PETSC_FALSE;
1232: #endif
1233:   PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_binary", &dump, NULL));
1234:   if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc));

1236:   if (viewer) {
1237:     for (PetscInt i = 1; i < n; i++) PetscCall(MatView(mglevels[i]->restrct, viewer));
1238:     for (PetscInt i = 0; i < n; i++) {
1239:       PetscCall(KSPGetPC(mglevels[i]->smoothd, &pc));
1240:       PetscCall(MatView(pc->mat, viewer));
1241:     }
1242:   }
1243:   PetscFunctionReturn(PETSC_SUCCESS);
1244: }

1246: PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels)
1247: {
1248:   PC_MG *mg = (PC_MG *)pc->data;

1250:   PetscFunctionBegin;
1251:   *levels = mg->nlevels;
1252:   PetscFunctionReturn(PETSC_SUCCESS);
1253: }

1255: /*@
1256:   PCMGGetLevels - Gets the number of levels to use with `PCMG`.

1258:   Not Collective

1260:   Input Parameter:
1261: . pc - the preconditioner context

1263:   Output Parameter:
1264: . levels - the number of levels

1266:   Level: advanced

1268: .seealso: [](ch_ksp), `PCMG`, `PCMGSetLevels()`
1269: @*/
1270: PetscErrorCode PCMGGetLevels(PC pc, PetscInt *levels)
1271: {
1272:   PetscFunctionBegin;
1274:   PetscAssertPointer(levels, 2);
1275:   *levels = 0;
1276:   PetscTryMethod(pc, "PCMGGetLevels_C", (PC, PetscInt *), (pc, levels));
1277:   PetscFunctionReturn(PETSC_SUCCESS);
1278: }

1280: /*@
1281:   PCMGGetGridComplexity - compute operator and grid complexity of the `PCMG` hierarchy

1283:   Input Parameter:
1284: . pc - the preconditioner context

1286:   Output Parameters:
1287: + gc - grid complexity = sum_i(n_i) / n_0
1288: - oc - operator complexity = sum_i(nnz_i) / nnz_0

1290:   Level: advanced

1292:   Note:
1293:   This is often call the operator complexity in multigrid literature

1295: .seealso: [](ch_ksp), `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`
1296: @*/
1297: PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc, PetscReal *oc)
1298: {
1299:   PC_MG         *mg       = (PC_MG *)pc->data;
1300:   PC_MG_Levels **mglevels = mg->levels;
1301:   PetscInt       N;
1302:   PetscLogDouble nnz0 = 0, sgc = 0, soc = 0, n0 = 0;
1303:   MatInfo        info;

1305:   PetscFunctionBegin;
1307:   if (gc) PetscAssertPointer(gc, 2);
1308:   if (oc) PetscAssertPointer(oc, 3);
1309:   if (!pc->setupcalled) {
1310:     if (gc) *gc = 0;
1311:     if (oc) *oc = 0;
1312:     PetscFunctionReturn(PETSC_SUCCESS);
1313:   }
1314:   PetscCheck(mg->nlevels > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MG has no levels");
1315:   for (PetscInt lev = 0; lev < mg->nlevels; lev++) {
1316:     Mat dB;
1317:     PetscCall(KSPGetOperators(mglevels[lev]->smoothd, NULL, &dB));
1318:     PetscCall(MatGetInfo(dB, MAT_GLOBAL_SUM, &info)); /* global reduction */
1319:     PetscCall(MatGetSize(dB, &N, NULL));
1320:     sgc += N;
1321:     soc += info.nz_used;
1322:     if (lev == mg->nlevels - 1) {
1323:       nnz0 = info.nz_used;
1324:       n0   = N;
1325:     }
1326:   }
1327:   PetscCheck(n0 > 0 && gc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number for grid points on finest level is not available");
1328:   *gc = (PetscReal)(sgc / n0);
1329:   if (nnz0 > 0 && oc) *oc = (PetscReal)(soc / nnz0);
1330:   PetscFunctionReturn(PETSC_SUCCESS);
1331: }

1333: /*@
1334:   PCMGSetType - Determines the form of multigrid to use, either
1335:   multiplicative, additive, full, or the Kaskade algorithm.

1337:   Logically Collective

1339:   Input Parameters:
1340: + pc   - the preconditioner context
1341: - form - multigrid form, one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE`

1343:   Options Database Key:
1344: . -pc_mg_type form - Sets form, one of multiplicative, additive, full, kaskade

1346:   Level: advanced

1348: .seealso: [](ch_ksp), `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGGetType()`, `PCMGCycleType`
1349: @*/
1350: PetscErrorCode PCMGSetType(PC pc, PCMGType form)
1351: {
1352:   PC_MG *mg = (PC_MG *)pc->data;

1354:   PetscFunctionBegin;
1357:   mg->am = form;
1358:   if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG;
1359:   else pc->ops->applyrichardson = NULL;
1360:   PetscFunctionReturn(PETSC_SUCCESS);
1361: }

1363: /*@
1364:   PCMGGetType - Finds the form of multigrid the `PCMG` is using  multiplicative, additive, full, or the Kaskade algorithm.

1366:   Logically Collective

1368:   Input Parameter:
1369: . pc - the preconditioner context

1371:   Output Parameter:
1372: . type - one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE`, `PCMGCycleType`

1374:   Level: advanced

1376: .seealso: [](ch_ksp), `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGSetType()`
1377: @*/
1378: PetscErrorCode PCMGGetType(PC pc, PCMGType *type)
1379: {
1380:   PC_MG *mg = (PC_MG *)pc->data;

1382:   PetscFunctionBegin;
1384:   *type = mg->am;
1385:   PetscFunctionReturn(PETSC_SUCCESS);
1386: }

1388: /*@
1389:   PCMGSetCycleType - Sets the type cycles to use.  Use `PCMGSetCycleTypeOnLevel()` for more
1390:   complicated cycling.

1392:   Logically Collective

1394:   Input Parameters:
1395: + pc - the multigrid context
1396: - n  - either `PC_MG_CYCLE_V` or `PC_MG_CYCLE_W`

1398:   Options Database Key:
1399: . -pc_mg_cycle_type (v|w) - provide the cycle desired

1401:   Level: advanced

1403: .seealso: [](ch_ksp), `PCMG`, `PCMGSetCycleTypeOnLevel()`, `PCMGType`, `PCMGCycleType`
1404: @*/
1405: PetscErrorCode PCMGSetCycleType(PC pc, PCMGCycleType n)
1406: {
1407:   PC_MG         *mg       = (PC_MG *)pc->data;
1408:   PC_MG_Levels **mglevels = mg->levels;
1409:   PetscInt       i, levels;

1411:   PetscFunctionBegin;
1414:   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling");
1415:   levels = mglevels[0]->levels;
1416:   for (i = 0; i < levels; i++) mglevels[i]->cycles = n;
1417:   PetscFunctionReturn(PETSC_SUCCESS);
1418: }

1420: /*@
1421:   PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step
1422:   of multigrid when `PCMGType` is `PC_MG_MULTIPLICATIVE`

1424:   Logically Collective

1426:   Input Parameters:
1427: + pc - the multigrid context
1428: - n  - number of cycles (default is 1)

1430:   Options Database Key:
1431: . -pc_mg_multiplicative_cycles n - set the number of cycles

1433:   Level: advanced

1435:   Note:
1436:   This is not associated with setting a v or w cycle, that is set with `PCMGSetCycleType()`

1438: .seealso: [](ch_ksp), `PCMGSetCycleTypeOnLevel()`, `PCMGSetCycleType()`, `PCMGCycleType`, `PCMGType`
1439: @*/
1440: PetscErrorCode PCMGMultiplicativeSetCycles(PC pc, PetscInt n)
1441: {
1442:   PC_MG *mg = (PC_MG *)pc->data;

1444:   PetscFunctionBegin;
1447:   mg->cyclesperpcapply = n;
1448:   PetscFunctionReturn(PETSC_SUCCESS);
1449: }

1451: /*
1452:    Since the finest level KSP shares the original matrix (of the entire system), it's preconditioner
1453:    should not be updated if the whole PC is supposed to reuse the preconditioner
1454: */
1455: static PetscErrorCode PCSetReusePreconditioner_MG(PC pc, PetscBool flag)
1456: {
1457:   PC_MG         *mg       = (PC_MG *)pc->data;
1458:   PC_MG_Levels **mglevels = mg->levels;
1459:   PetscInt       levels;
1460:   PC             tpc;

1462:   PetscFunctionBegin;
1465:   if (mglevels) {
1466:     levels = mglevels[0]->levels;
1467:     PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc));
1468:     tpc->reusepreconditioner = flag;
1469:     PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc));
1470:     tpc->reusepreconditioner = flag;
1471:   }
1472:   PetscFunctionReturn(PETSC_SUCCESS);
1473: }

1475: static PetscErrorCode PCMGSetGalerkin_MG(PC pc, PCMGGalerkinType use)
1476: {
1477:   PC_MG *mg = (PC_MG *)pc->data;

1479:   PetscFunctionBegin;
1480:   mg->galerkin = use;
1481:   PetscFunctionReturn(PETSC_SUCCESS);
1482: }

1484: /*@
1485:   PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the
1486:   finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i

1488:   Logically Collective

1490:   Input Parameters:
1491: + pc  - the multigrid context
1492: - use - one of `PC_MG_GALERKIN_BOTH`, `PC_MG_GALERKIN_PMAT`, `PC_MG_GALERKIN_MAT`, or `PC_MG_GALERKIN_NONE`

1494:   Options Database Key:
1495: . -pc_mg_galerkin (both|pmat|mat|none) - set the matrices to form via the Galerkin process

1497:   Level: intermediate

1499:   Note:
1500:   Some codes that use `PCMG` such as `PCGAMG` use Galerkin internally while constructing the hierarchy and thus do not
1501:   use the `PCMG` construction of the coarser grids.

1503: .seealso: [](ch_ksp), `PCMG`, `PCMGGetGalerkin()`, `PCMGGalerkinType`
1504: @*/
1505: PetscErrorCode PCMGSetGalerkin(PC pc, PCMGGalerkinType use)
1506: {
1507:   PetscFunctionBegin;
1509:   PetscTryMethod(pc, "PCMGSetGalerkin_C", (PC, PCMGGalerkinType), (pc, use));
1510:   PetscFunctionReturn(PETSC_SUCCESS);
1511: }

1513: /*@
1514:   PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. A_i-1 = r_i * A_i * p_i

1516:   Not Collective

1518:   Input Parameter:
1519: . pc - the multigrid context

1521:   Output Parameter:
1522: . galerkin - one of `PC_MG_GALERKIN_BOTH`,`PC_MG_GALERKIN_PMAT`,`PC_MG_GALERKIN_MAT`, `PC_MG_GALERKIN_NONE`, or `PC_MG_GALERKIN_EXTERNAL`

1524:   Level: intermediate

1526: .seealso: [](ch_ksp), `PCMG`, `PCMGSetGalerkin()`, `PCMGGalerkinType`
1527: @*/
1528: PetscErrorCode PCMGGetGalerkin(PC pc, PCMGGalerkinType *galerkin)
1529: {
1530:   PC_MG *mg = (PC_MG *)pc->data;

1532:   PetscFunctionBegin;
1534:   *galerkin = mg->galerkin;
1535:   PetscFunctionReturn(PETSC_SUCCESS);
1536: }

1538: static PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt)
1539: {
1540:   PC_MG *mg = (PC_MG *)pc->data;

1542:   PetscFunctionBegin;
1543:   mg->adaptInterpolation = adapt;
1544:   PetscFunctionReturn(PETSC_SUCCESS);
1545: }

1547: static PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt)
1548: {
1549:   PC_MG *mg = (PC_MG *)pc->data;

1551:   PetscFunctionBegin;
1552:   *adapt = mg->adaptInterpolation;
1553:   PetscFunctionReturn(PETSC_SUCCESS);
1554: }

1556: static PetscErrorCode PCMGSetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType ctype)
1557: {
1558:   PC_MG *mg = (PC_MG *)pc->data;

1560:   PetscFunctionBegin;
1561:   mg->adaptInterpolation = ctype != PCMG_ADAPT_NONE ? PETSC_TRUE : PETSC_FALSE;
1562:   mg->coarseSpaceType    = ctype;
1563:   PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_BOTH));
1564:   PetscFunctionReturn(PETSC_SUCCESS);
1565: }

1567: static PetscErrorCode PCMGGetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType *ctype)
1568: {
1569:   PC_MG *mg = (PC_MG *)pc->data;

1571:   PetscFunctionBegin;
1572:   *ctype = mg->coarseSpaceType;
1573:   PetscFunctionReturn(PETSC_SUCCESS);
1574: }

1576: static PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr)
1577: {
1578:   PC_MG *mg = (PC_MG *)pc->data;

1580:   PetscFunctionBegin;
1581:   mg->compatibleRelaxation = cr;
1582:   PetscFunctionReturn(PETSC_SUCCESS);
1583: }

1585: static PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr)
1586: {
1587:   PC_MG *mg = (PC_MG *)pc->data;

1589:   PetscFunctionBegin;
1590:   *cr = mg->compatibleRelaxation;
1591:   PetscFunctionReturn(PETSC_SUCCESS);
1592: }

1594: /*@
1595:   PCMGSetAdaptCoarseSpaceType - Set the type of adaptive coarse space.

1597:   Adapts or creates the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated.

1599:   Logically Collective

1601:   Input Parameters:
1602: + pc    - the multigrid context
1603: - ctype - the type of coarse space

1605:   Options Database Keys:
1606: + -pc_mg_adapt_interp_n nmodes          - The number of modes to use
1607: - -pc_mg_adapt_interp_coarse_space type - The type of coarse space: none, `polynomial`, `harmonic`, `eigenvector`, `generalized_eigenvector`, `gdsw`

1609:   Level: intermediate

1611:   Note:
1612:   Requires a `DM` with specific functionality be attached to the `PC`.

1614: .seealso: [](ch_ksp), `PCMG`, `PCMGCoarseSpaceType`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()`, `DM`
1615: @*/
1616: PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType ctype)
1617: {
1618:   PetscFunctionBegin;
1621:   PetscTryMethod(pc, "PCMGSetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType), (pc, ctype));
1622:   PetscFunctionReturn(PETSC_SUCCESS);
1623: }

1625: /*@
1626:   PCMGGetAdaptCoarseSpaceType - Get the type of adaptive coarse space.

1628:   Not Collective

1630:   Input Parameter:
1631: . pc - the multigrid context

1633:   Output Parameter:
1634: . ctype - the type of coarse space

1636:   Level: intermediate

1638: .seealso: [](ch_ksp), `PCMG`, `PCMGCoarseSpaceType`, `PCMGSetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()`
1639: @*/
1640: PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType *ctype)
1641: {
1642:   PetscFunctionBegin;
1644:   PetscAssertPointer(ctype, 2);
1645:   PetscUseMethod(pc, "PCMGGetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType *), (pc, ctype));
1646:   PetscFunctionReturn(PETSC_SUCCESS);
1647: }

1649: /* MATT: REMOVE? */
1650: /*@
1651:   PCMGSetAdaptInterpolation - Adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated.

1653:   Logically Collective

1655:   Input Parameters:
1656: + pc    - the multigrid context
1657: - adapt - flag for adaptation of the interpolator

1659:   Options Database Keys:
1660: + -pc_mg_adapt_interp                   - Turn on adaptation
1661: . -pc_mg_adapt_interp_n nmodes          - The number of modes to use, should be divisible by dimension
1662: - -pc_mg_adapt_interp_coarse_space type - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector

1664:   Level: intermediate

1666: .seealso: [](ch_ksp), `PCMG`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
1667: @*/
1668: PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt)
1669: {
1670:   PetscFunctionBegin;
1672:   PetscTryMethod(pc, "PCMGSetAdaptInterpolation_C", (PC, PetscBool), (pc, adapt));
1673:   PetscFunctionReturn(PETSC_SUCCESS);
1674: }

1676: /*@
1677:   PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh,
1678:   and thus accurately interpolated.

1680:   Not Collective

1682:   Input Parameter:
1683: . pc - the multigrid context

1685:   Output Parameter:
1686: . adapt - flag for adaptation of the interpolator

1688:   Level: intermediate

1690: .seealso: [](ch_ksp), `PCMG`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
1691: @*/
1692: PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt)
1693: {
1694:   PetscFunctionBegin;
1696:   PetscAssertPointer(adapt, 2);
1697:   PetscUseMethod(pc, "PCMGGetAdaptInterpolation_C", (PC, PetscBool *), (pc, adapt));
1698:   PetscFunctionReturn(PETSC_SUCCESS);
1699: }

1701: /*@
1702:   PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation.

1704:   Logically Collective

1706:   Input Parameters:
1707: + pc - the multigrid context
1708: - cr - flag for compatible relaxation

1710:   Options Database Key:
1711: . -pc_mg_adapt_cr - Turn on compatible relaxation

1713:   Level: intermediate

1715: .seealso: [](ch_ksp), `PCMG`, `PCMGGetAdaptCR()`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
1716: @*/
1717: PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr)
1718: {
1719:   PetscFunctionBegin;
1721:   PetscTryMethod(pc, "PCMGSetAdaptCR_C", (PC, PetscBool), (pc, cr));
1722:   PetscFunctionReturn(PETSC_SUCCESS);
1723: }

1725: /*@
1726:   PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation.

1728:   Not Collective

1730:   Input Parameter:
1731: . pc - the multigrid context

1733:   Output Parameter:
1734: . cr - flag for compatible relaxaion

1736:   Level: intermediate

1738: .seealso: [](ch_ksp), `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
1739: @*/
1740: PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr)
1741: {
1742:   PetscFunctionBegin;
1744:   PetscAssertPointer(cr, 2);
1745:   PetscUseMethod(pc, "PCMGGetAdaptCR_C", (PC, PetscBool *), (pc, cr));
1746:   PetscFunctionReturn(PETSC_SUCCESS);
1747: }

1749: /*@
1750:   PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use
1751:   on all levels.  Use `PCMGDistinctSmoothUp()` to create separate up and down smoothers if you want different numbers of
1752:   pre- and post-smoothing steps.

1754:   Logically Collective

1756:   Input Parameters:
1757: + pc - the multigrid context
1758: - n  - the number of smoothing steps

1760:   Options Database Key:
1761: . -mg_levels_ksp_max_it n - Sets number of pre and post-smoothing steps

1763:   Level: advanced

1765:   Note:
1766:   This does not set a value on the coarsest grid, since we assume that there is no separate smooth up on the coarsest grid.

1768: .seealso: [](ch_ksp), `PCMG`, `PCMGSetDistinctSmoothUp()`
1769: @*/
1770: PetscErrorCode PCMGSetNumberSmooth(PC pc, PetscInt n)
1771: {
1772:   PC_MG         *mg       = (PC_MG *)pc->data;
1773:   PC_MG_Levels **mglevels = mg->levels;
1774:   PetscInt       levels;

1776:   PetscFunctionBegin;
1779:   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling");
1780:   levels = mglevels[0]->levels;

1782:   for (PetscInt i = 1; i < levels; i++) {
1783:     PetscCall(KSPSetTolerances(mglevels[i]->smoothu, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, n));
1784:     PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, n));
1785:     mg->default_smoothu = n;
1786:     mg->default_smoothd = n;
1787:   }
1788:   PetscFunctionReturn(PETSC_SUCCESS);
1789: }

1791: /*@
1792:   PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate `KSP` from the down (pre) smoother on all levels
1793:   and adds the suffix _up to the options name

1795:   Logically Collective

1797:   Input Parameter:
1798: . pc - the preconditioner context

1800:   Options Database Key:
1801: . -pc_mg_distinct_smoothup (true|false) - use distinct smoothing objects

1803:   Level: advanced

1805:   Note:
1806:   This does not set a value on the coarsest grid, since we assume that there is no separate smooth up on the coarsest grid.

1808: .seealso: [](ch_ksp), `PCMG`, `PCMGSetNumberSmooth()`
1809: @*/
1810: PetscErrorCode PCMGSetDistinctSmoothUp(PC pc)
1811: {
1812:   PC_MG         *mg       = (PC_MG *)pc->data;
1813:   PC_MG_Levels **mglevels = mg->levels;
1814:   PetscInt       levels;
1815:   KSP            subksp;

1817:   PetscFunctionBegin;
1819:   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling");
1820:   levels = mglevels[0]->levels;

1822:   for (PetscInt i = 1; i < levels; i++) {
1823:     const char *prefix = NULL;
1824:     /* make sure smoother up and down are different */
1825:     PetscCall(PCMGGetSmootherUp(pc, i, &subksp));
1826:     PetscCall(KSPGetOptionsPrefix(mglevels[i]->smoothd, &prefix));
1827:     PetscCall(KSPSetOptionsPrefix(subksp, prefix));
1828:     PetscCall(KSPAppendOptionsPrefix(subksp, "up_"));
1829:   }
1830:   PetscFunctionReturn(PETSC_SUCCESS);
1831: }

1833: /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1834: static PetscErrorCode PCGetInterpolations_MG(PC pc, PetscInt *num_levels, Mat *interpolations[])
1835: {
1836:   PC_MG         *mg       = (PC_MG *)pc->data;
1837:   PC_MG_Levels **mglevels = mg->levels;
1838:   Mat           *mat;

1840:   PetscFunctionBegin;
1841:   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
1842:   PetscCall(PetscMalloc1(mg->nlevels, &mat));
1843:   for (PetscInt l = 1; l < mg->nlevels; l++) {
1844:     mat[l - 1] = mglevels[l]->interpolate;
1845:     PetscCall(PetscObjectReference((PetscObject)mat[l - 1]));
1846:   }
1847:   *num_levels     = mg->nlevels;
1848:   *interpolations = mat;
1849:   PetscFunctionReturn(PETSC_SUCCESS);
1850: }

1852: /* No new matrices are created, and the coarse operator matrices are the references to the original ones */
1853: static PetscErrorCode PCGetCoarseOperators_MG(PC pc, PetscInt *num_levels, Mat *coarseOperators[])
1854: {
1855:   PC_MG         *mg       = (PC_MG *)pc->data;
1856:   PC_MG_Levels **mglevels = mg->levels;
1857:   Mat           *mat;

1859:   PetscFunctionBegin;
1860:   PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling");
1861:   PetscCall(PetscMalloc1(mg->nlevels, &mat));
1862:   for (PetscInt l = 0; l < mg->nlevels - 1; l++) {
1863:     PetscCall(KSPGetOperators(mglevels[l]->smoothd, NULL, &mat[l]));
1864:     PetscCall(PetscObjectReference((PetscObject)mat[l]));
1865:   }
1866:   *num_levels      = mg->nlevels;
1867:   *coarseOperators = mat;
1868:   PetscFunctionReturn(PETSC_SUCCESS);
1869: }

1871: /*@C
1872:   PCMGRegisterCoarseSpaceConstructor -  Adds a method to the `PCMG` package for coarse space construction.

1874:   Not Collective, No Fortran Support

1876:   Input Parameters:
1877: + name     - name of the constructor
1878: - function - constructor routine, see `PCMGCoarseSpaceConstructorFn`

1880:   Level: advanced

1882:   Developer Notes:
1883:   This does not appear to be used anywhere

1885: .seealso: [](ch_ksp), `PCMGCoarseSpaceConstructorFn`, `PCMG`, `PCMGGetCoarseSpaceConstructor()`, `PCRegister()`
1886: @*/
1887: PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PCMGCoarseSpaceConstructorFn *function)
1888: {
1889:   PetscFunctionBegin;
1890:   PetscCall(PCInitializePackage());
1891:   PetscCall(PetscFunctionListAdd(&PCMGCoarseList, name, function));
1892:   PetscFunctionReturn(PETSC_SUCCESS);
1893: }

1895: /*@C
1896:   PCMGGetCoarseSpaceConstructor -  Returns the given coarse space construction method.

1898:   Not Collective, No Fortran Support

1900:   Input Parameter:
1901: . name - name of the constructor

1903:   Output Parameter:
1904: . function - constructor routine

1906:   Level: advanced

1908: .seealso: [](ch_ksp), `PCMGCoarseSpaceConstructorFn`, `PCMG`, `PCMGRegisterCoarseSpaceConstructor()`, `PCRegister()`
1909: @*/
1910: PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PCMGCoarseSpaceConstructorFn **function)
1911: {
1912:   PetscFunctionBegin;
1913:   PetscCall(PetscFunctionListFind(PCMGCoarseList, name, function));
1914:   PetscFunctionReturn(PETSC_SUCCESS);
1915: }

1917: /*MC
1918:    PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional information about the restriction/interpolation
1919:    operators using `PCMGSetInterpolation()` and/or `PCMGSetRestriction()`(and possibly the coarser grid matrices) or a `DM` that can provide such information.

1921:    Options Database Keys:
1922: +  -pc_mg_levels nlevels                              - number of levels including finest
1923: .  -pc_mg_cycle_type (v|w)                            - provide the cycle desired
1924: .  -pc_mg_type (additive|multiplicative|full|kaskade) - multiplicative is the default
1925: .  -pc_mg_log                                         - log information about time spent on each level of the solver
1926: .  -pc_mg_distinct_smoothup                           - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes)
1927: .  -pc_mg_galerkin (both|pmat|mat|none)               - use the Galerkin process to compute coarser operators, i.e., $A_{coarse} = R A_{fine} R^T$
1928: .  -pc_mg_multiplicative_cycles ncycles               - number of cycles to use as the preconditioner (defaults to 1)
1929: .  -pc_mg_dump_matlab                                 - dumps the matrices for each level and the restriction/interpolation matrices
1930:                                                         to a `PETSCVIEWERSOCKET` for reading from MATLAB.
1931: -  -pc_mg_dump_binary                                 - dumps the matrices for each level and the restriction/interpolation matrices
1932:                                                         to the binary output file called binaryoutput

1934:    Level: intermediate

1936:    Notes:
1937:    `PCMG` provides a general framework for implementing multigrid methods. Use `PCGAMG` for PETSc's algebraic multigrid preconditioner, `PCHYPRE` for hypre's
1938:    BoomerAMG algebraic multigrid, and `PCML` for Trilinos's ML preconditioner. `PCAMGX` provides access to NVIDIA's AmgX algebraic multigrid.

1940:    If you use `KSPSetDM()` (or `SNESSetDM()` or `TSSetDM()`) with an appropriate `DM`, such as `DMDA`, then `PCMG` will use the geometric information
1941:    from the `DM` to generate appropriate restriction and interpolation information and construct a geometric multigrid.

1943:    If you do not provide an appropriate `DM` and do not provide restriction or interpolation operators with `PCMGSetInterpolation()` and/or `PCMGSetRestriction()`,
1944:    then `PCMG` will run multigrid with only a single level (so not really multigrid).

1946:    The Krylov solver (if any) and preconditioner (smoother) and their parameters are controlled from the options database with the standard
1947:    options database keywords prefixed with `-mg_levels_` to affect all the levels but the coarsest, which is controlled with `-mg_coarse_`,
1948:    and the finest where `-mg_fine_` can override `-mg_levels_`.  One can set different preconditioners etc on specific levels with the prefix
1949:    `-mg_levels_n_` where `n` is the level number (zero being the coarse level. For example
1950: .vb
1951:    -mg_levels_ksp_type gmres -mg_levels_pc_type bjacobi -mg_coarse_pc_type svd -mg_levels_2_pc_type sor
1952: .ve
1953:    These options also work for controlling the smoothers etc inside `PCGAMG`

1955:    If one uses a Krylov method such `KSPGMRES` or `KSPCG` as the smoother than one must use `KSPFGMRES`, `KSPGCR`, or `KSPRICHARDSON` as the outer Krylov method

1957:    When run with a single level the smoother options are used on that level NOT the coarse grid solver options

1959:    When run with `KSPRICHARDSON` the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This
1960:    is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing
1961:    (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the
1962:    residual is computed at the end of each cycle.

1964: .seealso: [](sec_mg), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCEXOTIC`, `PCGAMG`, `PCML`, `PCHYPRE`,
1965:           `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `PCMGSetCycleType()`,
1966:           `PCMGSetDistinctSmoothUp()`, `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`,
1967:           `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`,
1968:           `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`,
1969:           `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()`
1970: M*/

1972: PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc)
1973: {
1974:   PC_MG *mg;

1976:   PetscFunctionBegin;
1977:   PetscCall(PetscNew(&mg));
1978:   pc->data               = mg;
1979:   mg->nlevels            = -1;
1980:   mg->am                 = PC_MG_MULTIPLICATIVE;
1981:   mg->galerkin           = PC_MG_GALERKIN_NONE;
1982:   mg->adaptInterpolation = PETSC_FALSE;
1983:   mg->Nc                 = -1;
1984:   mg->eigenvalue         = -1;

1986:   pc->useAmat = PETSC_TRUE;

1988:   pc->ops->apply             = PCApply_MG;
1989:   pc->ops->applytranspose    = PCApplyTranspose_MG;
1990:   pc->ops->matapply          = PCMatApply_MG;
1991:   pc->ops->matapplytranspose = PCMatApplyTranspose_MG;
1992:   pc->ops->setup             = PCSetUp_MG;
1993:   pc->ops->reset             = PCReset_MG;
1994:   pc->ops->destroy           = PCDestroy_MG;
1995:   pc->ops->setfromoptions    = PCSetFromOptions_MG;
1996:   pc->ops->view              = PCView_MG;

1998:   PetscCall(PetscObjectComposedDataRegister(&mg->eigenvalue));
1999:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", PCMGSetGalerkin_MG));
2000:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetReusePreconditioner_C", PCSetReusePreconditioner_MG));
2001:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", PCMGGetLevels_MG));
2002:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", PCMGSetLevels_MG));
2003:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", PCGetInterpolations_MG));
2004:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", PCGetCoarseOperators_MG));
2005:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", PCMGSetAdaptInterpolation_MG));
2006:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", PCMGGetAdaptInterpolation_MG));
2007:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", PCMGSetAdaptCR_MG));
2008:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", PCMGGetAdaptCR_MG));
2009:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", PCMGSetAdaptCoarseSpaceType_MG));
2010:   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", PCMGGetAdaptCoarseSpaceType_MG));
2011:   PetscFunctionReturn(PETSC_SUCCESS);
2012: }