1 /* 2 Defines the multigrid preconditioner interface. 3 */ 4 #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 5 #include <petsc/private/kspimpl.h> 6 #include <petscdm.h> 7 PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC, PetscBool *); 8 9 /* 10 Contains the list of registered coarse space construction routines 11 */ 12 PetscFunctionList PCMGCoarseList = NULL; 13 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; 19 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)); 47 48 /* if on finest level and have convergence criteria set */ 49 if (mglevels->level == mglevels->levels - 1 && mg->ttol && reason) { 50 PetscReal rnorm; 51 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 } 64 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) { 77 PetscCall(MatDuplicate(mgc->B, MAT_DO_NOT_COPY_VALUES, &mgc->X)); 78 } else { 79 PetscCall(MatZeroEntries(mgc->X)); 80 } 81 } else { 82 PetscCall(VecZeroEntries(mgc->x)); 83 } 84 while (cycles--) PetscCall(PCMGMCycle_Private(pc, mglevelsin - 1, transpose, matapp, reason)); 85 if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict, 0, 0, 0, 0)); 86 if (!transpose) { 87 if (matapp) PetscCall(MatMatInterpolateAdd(mglevels->interpolate, mgc->X, mglevels->X, &mglevels->X)); 88 else PetscCall(MatInterpolateAdd(mglevels->interpolate, mgc->x, mglevels->x, mglevels->x)); 89 } else { 90 PetscCall(MatInterpolateAdd(mglevels->restrct, mgc->x, mglevels->x, mglevels->x)); 91 } 92 if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict, 0, 0, 0, 0)); 93 if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve, 0, 0, 0, 0)); 94 if (!transpose) { 95 if (matapp) { 96 PetscCall(KSPMatSolve(mglevels->smoothu, mglevels->B, mglevels->X)); /* post smooth */ 97 PetscCall(KSPCheckSolve(mglevels->smoothu, pc, NULL)); 98 } else { 99 PetscCall(KSPSolve(mglevels->smoothu, mglevels->b, mglevels->x)); /* post smooth */ 100 PetscCall(KSPCheckSolve(mglevels->smoothu, pc, mglevels->x)); 101 } 102 } else { 103 PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported"); 104 PetscCall(KSPSolveTranspose(mglevels->smoothd, mglevels->b, mglevels->x)); /* post smooth */ 105 PetscCall(KSPCheckSolve(mglevels->smoothd, pc, mglevels->x)); 106 } 107 if (mglevels->cr) { 108 Mat crA; 109 110 PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported"); 111 /* TODO Turn on copy and turn off noisy if we have an exact solution 112 PetscCall(VecCopy(mglevels->x, mglevels->crx)); 113 PetscCall(VecCopy(mglevels->b, mglevels->crb)); */ 114 PetscCall(KSPGetOperators(mglevels->cr, &crA, NULL)); 115 PetscCall(KSPSetNoisy_Private(crA, mglevels->crx)); 116 PetscCall(KSPSolve(mglevels->cr, mglevels->crb, mglevels->crx)); /* compatible relaxation */ 117 PetscCall(KSPCheckSolve(mglevels->cr, pc, mglevels->crx)); 118 } 119 if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve, 0, 0, 0, 0)); 120 } 121 PetscFunctionReturn(PETSC_SUCCESS); 122 } 123 124 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) 125 { 126 PC_MG *mg = (PC_MG *)pc->data; 127 PC_MG_Levels **mglevels = mg->levels; 128 PC tpc; 129 PetscBool changeu, changed; 130 PetscInt levels = mglevels[0]->levels, i; 131 132 PetscFunctionBegin; 133 /* When the DM is supplying the matrix then it will not exist until here */ 134 for (i = 0; i < levels; i++) { 135 if (!mglevels[i]->A) { 136 PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL)); 137 PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A)); 138 } 139 } 140 141 PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc)); 142 PetscCall(PCPreSolveChangeRHS(tpc, &changed)); 143 PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc)); 144 PetscCall(PCPreSolveChangeRHS(tpc, &changeu)); 145 if (!changed && !changeu) { 146 PetscCall(VecDestroy(&mglevels[levels - 1]->b)); 147 mglevels[levels - 1]->b = b; 148 } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 149 if (!mglevels[levels - 1]->b) { 150 Vec *vec; 151 152 PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL)); 153 mglevels[levels - 1]->b = *vec; 154 PetscCall(PetscFree(vec)); 155 } 156 PetscCall(VecCopy(b, mglevels[levels - 1]->b)); 157 } 158 mglevels[levels - 1]->x = x; 159 160 mg->rtol = rtol; 161 mg->abstol = abstol; 162 mg->dtol = dtol; 163 if (rtol) { 164 /* compute initial residual norm for relative convergence test */ 165 PetscReal rnorm; 166 167 if (zeroguess) { 168 PetscCall(VecNorm(b, NORM_2, &rnorm)); 169 } else { 170 PetscCall((*mglevels[levels - 1]->residual)(mglevels[levels - 1]->A, b, x, w)); 171 PetscCall(VecNorm(w, NORM_2, &rnorm)); 172 } 173 mg->ttol = PetscMax(rtol * rnorm, abstol); 174 } else if (abstol) mg->ttol = abstol; 175 else mg->ttol = 0.0; 176 177 /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 178 stop prematurely due to small residual */ 179 for (i = 1; i < levels; i++) { 180 PetscCall(KSPSetTolerances(mglevels[i]->smoothu, 0, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 181 if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 182 /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 183 PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE)); 184 PetscCall(KSPSetTolerances(mglevels[i]->smoothd, 0, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT)); 185 } 186 } 187 188 *reason = PCRICHARDSON_NOT_SET; 189 for (i = 0; i < its; i++) { 190 PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, PETSC_FALSE, PETSC_FALSE, reason)); 191 if (*reason) break; 192 } 193 if (*reason == PCRICHARDSON_NOT_SET) *reason = PCRICHARDSON_CONVERGED_ITS; 194 *outits = i; 195 if (!changed && !changeu) mglevels[levels - 1]->b = NULL; 196 PetscFunctionReturn(PETSC_SUCCESS); 197 } 198 199 PetscErrorCode PCReset_MG(PC pc) 200 { 201 PC_MG *mg = (PC_MG *)pc->data; 202 PC_MG_Levels **mglevels = mg->levels; 203 PetscInt i, n; 204 205 PetscFunctionBegin; 206 if (mglevels) { 207 n = mglevels[0]->levels; 208 for (i = 0; i < n - 1; i++) { 209 PetscCall(VecDestroy(&mglevels[i + 1]->r)); 210 PetscCall(VecDestroy(&mglevels[i]->b)); 211 PetscCall(VecDestroy(&mglevels[i]->x)); 212 PetscCall(MatDestroy(&mglevels[i + 1]->R)); 213 PetscCall(MatDestroy(&mglevels[i]->B)); 214 PetscCall(MatDestroy(&mglevels[i]->X)); 215 PetscCall(VecDestroy(&mglevels[i]->crx)); 216 PetscCall(VecDestroy(&mglevels[i]->crb)); 217 PetscCall(MatDestroy(&mglevels[i + 1]->restrct)); 218 PetscCall(MatDestroy(&mglevels[i + 1]->interpolate)); 219 PetscCall(MatDestroy(&mglevels[i + 1]->inject)); 220 PetscCall(VecDestroy(&mglevels[i + 1]->rscale)); 221 } 222 PetscCall(VecDestroy(&mglevels[n - 1]->crx)); 223 PetscCall(VecDestroy(&mglevels[n - 1]->crb)); 224 /* this is not null only if the smoother on the finest level 225 changes the rhs during PreSolve */ 226 PetscCall(VecDestroy(&mglevels[n - 1]->b)); 227 PetscCall(MatDestroy(&mglevels[n - 1]->B)); 228 229 for (i = 0; i < n; i++) { 230 PetscCall(MatDestroy(&mglevels[i]->coarseSpace)); 231 PetscCall(MatDestroy(&mglevels[i]->A)); 232 if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPReset(mglevels[i]->smoothd)); 233 PetscCall(KSPReset(mglevels[i]->smoothu)); 234 if (mglevels[i]->cr) PetscCall(KSPReset(mglevels[i]->cr)); 235 } 236 mg->Nc = 0; 237 } 238 PetscFunctionReturn(PETSC_SUCCESS); 239 } 240 241 /* Implementing CR 242 243 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 245 Inj^T (Inj Inj^T)^{-1} Inj 246 247 and if Inj is a VecScatter, as it is now in PETSc, we have 248 249 Inj^T Inj 250 251 and 252 253 S = I - Inj^T Inj 254 255 since 256 257 Inj S = Inj - (Inj Inj^T) Inj = 0. 258 259 Brannick suggests 260 261 A \to S^T A S \qquad\mathrm{and}\qquad M \to S^T M S 262 263 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 265 M^{-1} A \to S M^{-1} A S 266 267 In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left. 268 269 Check: || Inj P - I ||_F < tol 270 Check: In general, Inj Inj^T = I 271 */ 272 273 typedef struct { 274 PC mg; /* The PCMG object */ 275 PetscInt l; /* The multigrid level for this solver */ 276 Mat Inj; /* The injection matrix */ 277 Mat S; /* I - Inj^T Inj */ 278 } CRContext; 279 280 static PetscErrorCode CRSetup_Private(PC pc) 281 { 282 CRContext *ctx; 283 Mat It; 284 285 PetscFunctionBeginUser; 286 PetscCall(PCShellGetContext(pc, &ctx)); 287 PetscCall(PCMGGetInjection(ctx->mg, ctx->l, &It)); 288 PetscCheck(It, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG"); 289 PetscCall(MatCreateTranspose(It, &ctx->Inj)); 290 PetscCall(MatCreateNormal(ctx->Inj, &ctx->S)); 291 PetscCall(MatScale(ctx->S, -1.0)); 292 PetscCall(MatShift(ctx->S, 1.0)); 293 PetscFunctionReturn(PETSC_SUCCESS); 294 } 295 296 static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y) 297 { 298 CRContext *ctx; 299 300 PetscFunctionBeginUser; 301 PetscCall(PCShellGetContext(pc, &ctx)); 302 PetscCall(MatMult(ctx->S, x, y)); 303 PetscFunctionReturn(PETSC_SUCCESS); 304 } 305 306 static PetscErrorCode CRDestroy_Private(PC pc) 307 { 308 CRContext *ctx; 309 310 PetscFunctionBeginUser; 311 PetscCall(PCShellGetContext(pc, &ctx)); 312 PetscCall(MatDestroy(&ctx->Inj)); 313 PetscCall(MatDestroy(&ctx->S)); 314 PetscCall(PetscFree(ctx)); 315 PetscCall(PCShellSetContext(pc, NULL)); 316 PetscFunctionReturn(PETSC_SUCCESS); 317 } 318 319 static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr) 320 { 321 CRContext *ctx; 322 323 PetscFunctionBeginUser; 324 PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), cr)); 325 PetscCall(PetscObjectSetName((PetscObject)*cr, "S (complementary projector to injection)")); 326 PetscCall(PetscCalloc1(1, &ctx)); 327 ctx->mg = pc; 328 ctx->l = l; 329 PetscCall(PCSetType(*cr, PCSHELL)); 330 PetscCall(PCShellSetContext(*cr, ctx)); 331 PetscCall(PCShellSetApply(*cr, CRApply_Private)); 332 PetscCall(PCShellSetSetUp(*cr, CRSetup_Private)); 333 PetscCall(PCShellSetDestroy(*cr, CRDestroy_Private)); 334 PetscFunctionReturn(PETSC_SUCCESS); 335 } 336 337 PETSC_EXTERN PetscErrorCode PetscOptionsFindPairPrefix_Private(PetscOptions, const char[], const char[], const char *[], const char *[], PetscBool *); 338 339 PetscErrorCode PCMGSetLevels_MG(PC pc, PetscInt levels, MPI_Comm *comms) 340 { 341 PC_MG *mg = (PC_MG *)pc->data; 342 MPI_Comm comm; 343 PC_MG_Levels **mglevels = mg->levels; 344 PCMGType mgtype = mg->am; 345 PetscInt mgctype = (PetscInt)PC_MG_CYCLE_V; 346 PetscInt i; 347 PetscMPIInt size; 348 const char *prefix; 349 PC ipc; 350 PetscInt n; 351 352 PetscFunctionBegin; 353 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 354 PetscValidLogicalCollectiveInt(pc, levels, 2); 355 if (mg->nlevels == levels) PetscFunctionReturn(PETSC_SUCCESS); 356 PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 357 if (mglevels) { 358 mgctype = mglevels[0]->cycles; 359 /* changing the number of levels so free up the previous stuff */ 360 PetscCall(PCReset_MG(pc)); 361 n = mglevels[0]->levels; 362 for (i = 0; i < n; i++) { 363 if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPDestroy(&mglevels[i]->smoothd)); 364 PetscCall(KSPDestroy(&mglevels[i]->smoothu)); 365 PetscCall(KSPDestroy(&mglevels[i]->cr)); 366 PetscCall(PetscFree(mglevels[i])); 367 } 368 PetscCall(PetscFree(mg->levels)); 369 } 370 371 mg->nlevels = levels; 372 373 PetscCall(PetscMalloc1(levels, &mglevels)); 374 375 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 376 377 mg->stageApply = 0; 378 for (i = 0; i < levels; i++) { 379 PetscCall(PetscNew(&mglevels[i])); 380 381 mglevels[i]->level = i; 382 mglevels[i]->levels = levels; 383 mglevels[i]->cycles = mgctype; 384 mg->default_smoothu = 2; 385 mg->default_smoothd = 2; 386 mglevels[i]->eventsmoothsetup = 0; 387 mglevels[i]->eventsmoothsolve = 0; 388 mglevels[i]->eventresidual = 0; 389 mglevels[i]->eventinterprestrict = 0; 390 391 if (comms) comm = comms[i]; 392 if (comm != MPI_COMM_NULL) { 393 PetscCall(KSPCreate(comm, &mglevels[i]->smoothd)); 394 PetscCall(KSPSetNestLevel(mglevels[i]->smoothd, pc->kspnestlevel)); 395 PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->smoothd, pc->erroriffailure)); 396 PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd, (PetscObject)pc, levels - i)); 397 PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd, prefix)); 398 PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level)); 399 if (i == 0 && levels > 1) { // coarse grid 400 PetscCall(KSPAppendOptionsPrefix(mglevels[0]->smoothd, "mg_coarse_")); 401 402 /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 403 PetscCall(KSPSetType(mglevels[0]->smoothd, KSPPREONLY)); 404 PetscCall(KSPGetPC(mglevels[0]->smoothd, &ipc)); 405 PetscCallMPI(MPI_Comm_size(comm, &size)); 406 if (size > 1) { 407 PetscCall(PCSetType(ipc, PCREDUNDANT)); 408 } else { 409 PetscCall(PCSetType(ipc, PCLU)); 410 } 411 PetscCall(PCFactorSetShiftType(ipc, MAT_SHIFT_INBLOCKS)); 412 } else { 413 char tprefix[128]; 414 415 PetscCall(KSPSetType(mglevels[i]->smoothd, KSPCHEBYSHEV)); 416 PetscCall(KSPSetConvergenceTest(mglevels[i]->smoothd, KSPConvergedSkip, NULL, NULL)); 417 PetscCall(KSPSetNormType(mglevels[i]->smoothd, KSP_NORM_NONE)); 418 PetscCall(KSPGetPC(mglevels[i]->smoothd, &ipc)); 419 PetscCall(PCSetType(ipc, PCSOR)); 420 PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, mg->default_smoothd)); 421 422 if (i == levels - 1 && levels > 1) { // replace 'mg_finegrid_' with 'mg_levels_X_' 423 PetscBool set; 424 425 PetscCall(PetscOptionsFindPairPrefix_Private(((PetscObject)mglevels[i]->smoothd)->options, ((PetscObject)mglevels[i]->smoothd)->prefix, "-mg_fine_", NULL, NULL, &set)); 426 if (set) { 427 if (prefix) PetscCall(PetscSNPrintf(tprefix, 128, "%smg_fine_", prefix)); 428 else PetscCall(PetscSNPrintf(tprefix, 128, "mg_fine_")); 429 PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd, tprefix)); 430 } else { 431 PetscCall(PetscSNPrintf(tprefix, 128, "mg_levels_%" PetscInt_FMT "_", i)); 432 PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd, tprefix)); 433 } 434 } else { 435 PetscCall(PetscSNPrintf(tprefix, 128, "mg_levels_%" PetscInt_FMT "_", i)); 436 PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd, tprefix)); 437 } 438 } 439 } 440 mglevels[i]->smoothu = mglevels[i]->smoothd; 441 mg->rtol = 0.0; 442 mg->abstol = 0.0; 443 mg->dtol = 0.0; 444 mg->ttol = 0.0; 445 mg->cyclesperpcapply = 1; 446 } 447 mg->levels = mglevels; 448 PetscCall(PCMGSetType(pc, mgtype)); 449 PetscFunctionReturn(PETSC_SUCCESS); 450 } 451 452 /*@C 453 PCMGSetLevels - Sets the number of levels to use with `PCMG`. 454 Must be called before any other `PCMG` routine. 455 456 Logically Collective 457 458 Input Parameters: 459 + pc - the preconditioner context 460 . levels - the number of levels 461 - comms - optional communicators for each level; this is to allow solving the coarser problems 462 on smaller sets of processes. For processes that are not included in the computation 463 you must pass `MPI_COMM_NULL`. Use comms = `NULL` to specify that all processes 464 should participate in each level of problem. 465 466 Level: intermediate 467 468 Notes: 469 If the number of levels is one then the multigrid uses the `-mg_levels` prefix 470 for setting the level options rather than the `-mg_coarse` or `-mg_fine` prefix. 471 472 You can free the information in comms after this routine is called. 473 474 The array of MPI communicators must contain `MPI_COMM_NULL` for those ranks that at each level 475 are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on 476 the two levels, rank 0 in the original communicator will pass in an array of 2 communicators 477 of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators 478 the first of size 2 and the second of value `MPI_COMM_NULL` since the rank 1 does not participate 479 in the coarse grid solve. 480 481 Since each coarser level may have a new `MPI_Comm` with fewer ranks than the previous, one 482 must take special care in providing the restriction and interpolation operation. We recommend 483 providing these as two step operations; first perform a standard restriction or interpolation on 484 the full number of ranks for that level and then use an MPI call to copy the resulting vector 485 array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, note in both 486 cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and 487 receives or `MPI_AlltoAllv()` could be used to do the reshuffling of the vector entries. 488 489 Fortran Notes: 490 Use comms = `PETSC_NULL_MPI_COMM` as the equivalent of `NULL` in the C interface. Note `PETSC_NULL_MPI_COMM` 491 is not `MPI_COMM_NULL`. It is more like `PETSC_NULL_INTEGER`, `PETSC_NULL_REAL` etc. 492 493 .seealso: [](ch_ksp), `PCMGSetType()`, `PCMGGetLevels()` 494 @*/ 495 PetscErrorCode PCMGSetLevels(PC pc, PetscInt levels, MPI_Comm *comms) 496 { 497 PetscFunctionBegin; 498 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 499 if (comms) PetscAssertPointer(comms, 3); 500 PetscTryMethod(pc, "PCMGSetLevels_C", (PC, PetscInt, MPI_Comm *), (pc, levels, comms)); 501 PetscFunctionReturn(PETSC_SUCCESS); 502 } 503 504 PetscErrorCode PCDestroy_MG(PC pc) 505 { 506 PC_MG *mg = (PC_MG *)pc->data; 507 PC_MG_Levels **mglevels = mg->levels; 508 PetscInt i, n; 509 510 PetscFunctionBegin; 511 PetscCall(PCReset_MG(pc)); 512 if (mglevels) { 513 n = mglevels[0]->levels; 514 for (i = 0; i < n; i++) { 515 if (mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPDestroy(&mglevels[i]->smoothd)); 516 PetscCall(KSPDestroy(&mglevels[i]->smoothu)); 517 PetscCall(KSPDestroy(&mglevels[i]->cr)); 518 PetscCall(PetscFree(mglevels[i])); 519 } 520 PetscCall(PetscFree(mg->levels)); 521 } 522 PetscCall(PetscFree(pc->data)); 523 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL)); 524 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL)); 525 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", NULL)); 526 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetReusePreconditioner_C", NULL)); 527 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", NULL)); 528 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", NULL)); 529 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", NULL)); 530 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", NULL)); 531 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", NULL)); 532 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", NULL)); 533 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", NULL)); 534 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", NULL)); 535 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", NULL)); 536 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", NULL)); 537 PetscFunctionReturn(PETSC_SUCCESS); 538 } 539 540 /* 541 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 542 or full cycle of multigrid. 543 544 Note: 545 A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 546 */ 547 static PetscErrorCode PCApply_MG_Internal(PC pc, Vec b, Vec x, Mat B, Mat X, PetscBool transpose) 548 { 549 PC_MG *mg = (PC_MG *)pc->data; 550 PC_MG_Levels **mglevels = mg->levels; 551 PC tpc; 552 PetscInt levels = mglevels[0]->levels, i; 553 PetscBool changeu, changed, matapp; 554 555 PetscFunctionBegin; 556 matapp = (PetscBool)(B && X); 557 if (mg->stageApply) PetscCall(PetscLogStagePush(mg->stageApply)); 558 /* When the DM is supplying the matrix then it will not exist until here */ 559 for (i = 0; i < levels; i++) { 560 if (!mglevels[i]->A) { 561 PetscCall(KSPGetOperators(mglevels[i]->smoothu, &mglevels[i]->A, NULL)); 562 PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A)); 563 } 564 } 565 566 PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc)); 567 PetscCall(PCPreSolveChangeRHS(tpc, &changed)); 568 PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc)); 569 PetscCall(PCPreSolveChangeRHS(tpc, &changeu)); 570 if (!changeu && !changed) { 571 if (matapp) { 572 PetscCall(MatDestroy(&mglevels[levels - 1]->B)); 573 mglevels[levels - 1]->B = B; 574 } else { 575 PetscCall(VecDestroy(&mglevels[levels - 1]->b)); 576 mglevels[levels - 1]->b = b; 577 } 578 } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 579 if (matapp) { 580 if (mglevels[levels - 1]->B) { 581 PetscInt N1, N2; 582 PetscBool flg; 583 584 PetscCall(MatGetSize(mglevels[levels - 1]->B, NULL, &N1)); 585 PetscCall(MatGetSize(B, NULL, &N2)); 586 PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels - 1]->B, ((PetscObject)B)->type_name, &flg)); 587 if (N1 != N2 || !flg) PetscCall(MatDestroy(&mglevels[levels - 1]->B)); 588 } 589 if (!mglevels[levels - 1]->B) { 590 PetscCall(MatDuplicate(B, MAT_COPY_VALUES, &mglevels[levels - 1]->B)); 591 } else { 592 PetscCall(MatCopy(B, mglevels[levels - 1]->B, SAME_NONZERO_PATTERN)); 593 } 594 } else { 595 if (!mglevels[levels - 1]->b) { 596 Vec *vec; 597 598 PetscCall(KSPCreateVecs(mglevels[levels - 1]->smoothd, 1, &vec, 0, NULL)); 599 mglevels[levels - 1]->b = *vec; 600 PetscCall(PetscFree(vec)); 601 } 602 PetscCall(VecCopy(b, mglevels[levels - 1]->b)); 603 } 604 } 605 if (matapp) { 606 mglevels[levels - 1]->X = X; 607 } else { 608 mglevels[levels - 1]->x = x; 609 } 610 611 /* If coarser Xs are present, it means we have already block applied the PC at least once 612 Reset operators if sizes/type do no match */ 613 if (matapp && levels > 1 && mglevels[levels - 2]->X) { 614 PetscInt Xc, Bc; 615 PetscBool flg; 616 617 PetscCall(MatGetSize(mglevels[levels - 2]->X, NULL, &Xc)); 618 PetscCall(MatGetSize(mglevels[levels - 1]->B, NULL, &Bc)); 619 PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels - 2]->X, ((PetscObject)mglevels[levels - 1]->X)->type_name, &flg)); 620 if (Xc != Bc || !flg) { 621 PetscCall(MatDestroy(&mglevels[levels - 1]->R)); 622 for (i = 0; i < levels - 1; i++) { 623 PetscCall(MatDestroy(&mglevels[i]->R)); 624 PetscCall(MatDestroy(&mglevels[i]->B)); 625 PetscCall(MatDestroy(&mglevels[i]->X)); 626 } 627 } 628 } 629 630 if (mg->am == PC_MG_MULTIPLICATIVE) { 631 if (matapp) PetscCall(MatZeroEntries(X)); 632 else PetscCall(VecZeroEntries(x)); 633 for (i = 0; i < mg->cyclesperpcapply; i++) PetscCall(PCMGMCycle_Private(pc, mglevels + levels - 1, transpose, matapp, NULL)); 634 } else if (mg->am == PC_MG_ADDITIVE) { 635 PetscCall(PCMGACycle_Private(pc, mglevels, transpose, matapp)); 636 } else if (mg->am == PC_MG_KASKADE) { 637 PetscCall(PCMGKCycle_Private(pc, mglevels, transpose, matapp)); 638 } else { 639 PetscCall(PCMGFCycle_Private(pc, mglevels, transpose, matapp)); 640 } 641 if (mg->stageApply) PetscCall(PetscLogStagePop()); 642 if (!changeu && !changed) { 643 if (matapp) { 644 mglevels[levels - 1]->B = NULL; 645 } else { 646 mglevels[levels - 1]->b = NULL; 647 } 648 } 649 PetscFunctionReturn(PETSC_SUCCESS); 650 } 651 652 static PetscErrorCode PCApply_MG(PC pc, Vec b, Vec x) 653 { 654 PetscFunctionBegin; 655 PetscCall(PCApply_MG_Internal(pc, b, x, NULL, NULL, PETSC_FALSE)); 656 PetscFunctionReturn(PETSC_SUCCESS); 657 } 658 659 static PetscErrorCode PCApplyTranspose_MG(PC pc, Vec b, Vec x) 660 { 661 PetscFunctionBegin; 662 PetscCall(PCApply_MG_Internal(pc, b, x, NULL, NULL, PETSC_TRUE)); 663 PetscFunctionReturn(PETSC_SUCCESS); 664 } 665 666 static PetscErrorCode PCMatApply_MG(PC pc, Mat b, Mat x) 667 { 668 PetscFunctionBegin; 669 PetscCall(PCApply_MG_Internal(pc, NULL, NULL, b, x, PETSC_FALSE)); 670 PetscFunctionReturn(PETSC_SUCCESS); 671 } 672 673 static PetscErrorCode PCMatApplyTranspose_MG(PC pc, Mat b, Mat x) 674 { 675 PetscFunctionBegin; 676 PetscCall(PCApply_MG_Internal(pc, NULL, NULL, b, x, PETSC_TRUE)); 677 PetscFunctionReturn(PETSC_SUCCESS); 678 } 679 680 PetscErrorCode PCSetFromOptions_MG(PC pc, PetscOptionItems PetscOptionsObject) 681 { 682 PetscInt levels, cycles; 683 PetscBool flg, flg2; 684 PC_MG *mg = (PC_MG *)pc->data; 685 PC_MG_Levels **mglevels; 686 PCMGType mgtype; 687 PCMGCycleType mgctype; 688 PCMGGalerkinType gtype; 689 PCMGCoarseSpaceType coarseSpaceType; 690 691 PetscFunctionBegin; 692 levels = PetscMax(mg->nlevels, 1); 693 PetscOptionsHeadBegin(PetscOptionsObject, "Multigrid options"); 694 PetscCall(PetscOptionsInt("-pc_mg_levels", "Number of Levels", "PCMGSetLevels", levels, &levels, &flg)); 695 if (!flg && !mg->levels && pc->dm) { 696 PetscCall(DMGetRefineLevel(pc->dm, &levels)); 697 levels++; 698 mg->usedmfornumberoflevels = PETSC_TRUE; 699 } 700 PetscCall(PCMGSetLevels(pc, levels, NULL)); 701 mglevels = mg->levels; 702 703 mgctype = (PCMGCycleType)mglevels[0]->cycles; 704 PetscCall(PetscOptionsEnum("-pc_mg_cycle_type", "V cycle or for W-cycle", "PCMGSetCycleType", PCMGCycleTypes, (PetscEnum)mgctype, (PetscEnum *)&mgctype, &flg)); 705 if (flg) PetscCall(PCMGSetCycleType(pc, mgctype)); 706 coarseSpaceType = mg->coarseSpaceType; 707 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)); 708 if (flg) PetscCall(PCMGSetAdaptCoarseSpaceType(pc, coarseSpaceType)); 709 PetscCall(PetscOptionsInt("-pc_mg_adapt_interp_n", "Size of the coarse space for adaptive interpolation", "PCMGSetCoarseSpace", mg->Nc, &mg->Nc, &flg)); 710 PetscCall(PetscOptionsBool("-pc_mg_mesp_monitor", "Monitor the multilevel eigensolver", "PCMGSetAdaptInterpolation", PETSC_FALSE, &mg->mespMonitor, &flg)); 711 flg2 = PETSC_FALSE; 712 PetscCall(PetscOptionsBool("-pc_mg_adapt_cr", "Monitor coarse space quality using Compatible Relaxation (CR)", "PCMGSetAdaptCR", PETSC_FALSE, &flg2, &flg)); 713 if (flg) PetscCall(PCMGSetAdaptCR(pc, flg2)); 714 flg = PETSC_FALSE; 715 PetscCall(PetscOptionsBool("-pc_mg_distinct_smoothup", "Create separate smoothup KSP and append the prefix _up", "PCMGSetDistinctSmoothUp", PETSC_FALSE, &flg, NULL)); 716 if (flg) PetscCall(PCMGSetDistinctSmoothUp(pc)); 717 PetscCall(PetscOptionsEnum("-pc_mg_galerkin", "Use Galerkin process to compute coarser operators", "PCMGSetGalerkin", PCMGGalerkinTypes, (PetscEnum)mg->galerkin, (PetscEnum *)>ype, &flg)); 718 if (flg) PetscCall(PCMGSetGalerkin(pc, gtype)); 719 mgtype = mg->am; 720 PetscCall(PetscOptionsEnum("-pc_mg_type", "Multigrid type", "PCMGSetType", PCMGTypes, (PetscEnum)mgtype, (PetscEnum *)&mgtype, &flg)); 721 if (flg) PetscCall(PCMGSetType(pc, mgtype)); 722 if (mg->am == PC_MG_MULTIPLICATIVE) { 723 PetscCall(PetscOptionsInt("-pc_mg_multiplicative_cycles", "Number of cycles for each preconditioner step", "PCMGMultiplicativeSetCycles", mg->cyclesperpcapply, &cycles, &flg)); 724 if (flg) PetscCall(PCMGMultiplicativeSetCycles(pc, cycles)); 725 } 726 flg = PETSC_FALSE; 727 PetscCall(PetscOptionsBool("-pc_mg_log", "Log times for each multigrid level", "None", flg, &flg, NULL)); 728 if (flg) { 729 PetscInt i; 730 char eventname[128]; 731 732 levels = mglevels[0]->levels; 733 for (i = 0; i < levels; i++) { 734 PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGSetup Level %" PetscInt_FMT, i)); 735 PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventsmoothsetup)); 736 PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGSmooth Level %" PetscInt_FMT, i)); 737 PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventsmoothsolve)); 738 if (i) { 739 PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGResid Level %" PetscInt_FMT, i)); 740 PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventresidual)); 741 PetscCall(PetscSNPrintf(eventname, PETSC_STATIC_ARRAY_LENGTH(eventname), "MGInterp Level %" PetscInt_FMT, i)); 742 PetscCall(PetscLogEventRegister(eventname, ((PetscObject)pc)->classid, &mglevels[i]->eventinterprestrict)); 743 } 744 } 745 746 if (PetscDefined(USE_LOG)) { 747 const char sname[] = "MG Apply"; 748 749 PetscCall(PetscLogStageGetId(sname, &mg->stageApply)); 750 if (mg->stageApply < 0) PetscCall(PetscLogStageRegister(sname, &mg->stageApply)); 751 } 752 } 753 PetscOptionsHeadEnd(); 754 /* Check option consistency */ 755 PetscCall(PCMGGetGalerkin(pc, >ype)); 756 PetscCall(PCMGGetAdaptInterpolation(pc, &flg)); 757 PetscCheck(!flg || !(gtype >= PC_MG_GALERKIN_NONE), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator"); 758 PetscFunctionReturn(PETSC_SUCCESS); 759 } 760 761 const char *const PCMGTypes[] = {"MULTIPLICATIVE", "ADDITIVE", "FULL", "KASKADE", "PCMGType", "PC_MG", NULL}; 762 const char *const PCMGCycleTypes[] = {"invalid", "v", "w", "PCMGCycleType", "PC_MG_CYCLE", NULL}; 763 const char *const PCMGGalerkinTypes[] = {"both", "pmat", "mat", "none", "external", "PCMGGalerkinType", "PC_MG_GALERKIN", NULL}; 764 const char *const PCMGCoarseSpaceTypes[] = {"none", "polynomial", "harmonic", "eigenvector", "generalized_eigenvector", "gdsw", "PCMGCoarseSpaceType", "PCMG_ADAPT_NONE", NULL}; 765 766 #include <petscdraw.h> 767 PetscErrorCode PCView_MG(PC pc, PetscViewer viewer) 768 { 769 PC_MG *mg = (PC_MG *)pc->data; 770 PC_MG_Levels **mglevels = mg->levels; 771 PetscInt levels = mglevels ? mglevels[0]->levels : 0, i; 772 PetscBool isascii, isbinary, isdraw; 773 774 PetscFunctionBegin; 775 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 776 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 777 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 778 if (isascii) { 779 const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 780 PetscCall(PetscViewerASCIIPrintf(viewer, " type is %s, levels=%" PetscInt_FMT " cycles=%s\n", PCMGTypes[mg->am], levels, cyclename)); 781 if (mg->am == PC_MG_MULTIPLICATIVE) PetscCall(PetscViewerASCIIPrintf(viewer, " Cycles per PCApply=%" PetscInt_FMT "\n", mg->cyclesperpcapply)); 782 if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 783 PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid matrices\n")); 784 } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 785 PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid matrices for pmat\n")); 786 } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 787 PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid matrices for mat\n")); 788 } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 789 PetscCall(PetscViewerASCIIPrintf(viewer, " Using externally compute Galerkin coarse grid matrices\n")); 790 } else { 791 PetscCall(PetscViewerASCIIPrintf(viewer, " Not using Galerkin computed coarse grid matrices\n")); 792 } 793 if (mg->view) PetscCall((*mg->view)(pc, viewer)); 794 for (i = 0; i < levels; i++) { 795 if (i) { 796 PetscCall(PetscViewerASCIIPrintf(viewer, "Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n", i)); 797 } else { 798 PetscCall(PetscViewerASCIIPrintf(viewer, "Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n", i)); 799 } 800 PetscCall(PetscViewerASCIIPushTab(viewer)); 801 PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 802 PetscCall(PetscViewerASCIIPopTab(viewer)); 803 if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 804 PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) same as down solver (pre-smoother)\n")); 805 } else if (i) { 806 PetscCall(PetscViewerASCIIPrintf(viewer, "Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n", i)); 807 PetscCall(PetscViewerASCIIPushTab(viewer)); 808 PetscCall(KSPView(mglevels[i]->smoothu, viewer)); 809 PetscCall(PetscViewerASCIIPopTab(viewer)); 810 } 811 if (i && mglevels[i]->cr) { 812 PetscCall(PetscViewerASCIIPrintf(viewer, "CR solver on level %" PetscInt_FMT " -------------------------------\n", i)); 813 PetscCall(PetscViewerASCIIPushTab(viewer)); 814 PetscCall(KSPView(mglevels[i]->cr, viewer)); 815 PetscCall(PetscViewerASCIIPopTab(viewer)); 816 } 817 } 818 } else if (isbinary) { 819 for (i = levels - 1; i >= 0; i--) { 820 PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 821 if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) PetscCall(KSPView(mglevels[i]->smoothu, viewer)); 822 } 823 } else if (isdraw) { 824 PetscDraw draw; 825 PetscReal x, w, y, bottom, th; 826 PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 827 PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 828 PetscCall(PetscDrawStringGetSize(draw, NULL, &th)); 829 bottom = y - th; 830 for (i = levels - 1; i >= 0; i--) { 831 if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 832 PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 833 PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 834 PetscCall(PetscDrawPopCurrentPoint(draw)); 835 } else { 836 w = 0.5 * PetscMin(1.0 - x, x); 837 PetscCall(PetscDrawPushCurrentPoint(draw, x + w, bottom)); 838 PetscCall(KSPView(mglevels[i]->smoothd, viewer)); 839 PetscCall(PetscDrawPopCurrentPoint(draw)); 840 PetscCall(PetscDrawPushCurrentPoint(draw, x - w, bottom)); 841 PetscCall(KSPView(mglevels[i]->smoothu, viewer)); 842 PetscCall(PetscDrawPopCurrentPoint(draw)); 843 } 844 PetscCall(PetscDrawGetBoundingBox(draw, NULL, &bottom, NULL, NULL)); 845 bottom -= th; 846 } 847 } 848 PetscFunctionReturn(PETSC_SUCCESS); 849 } 850 851 #include <petsc/private/kspimpl.h> 852 853 /* 854 Calls setup for the KSP on each level 855 */ 856 PetscErrorCode PCSetUp_MG(PC pc) 857 { 858 PC_MG *mg = (PC_MG *)pc->data; 859 PC_MG_Levels **mglevels = mg->levels; 860 PetscInt i, n; 861 PC cpc; 862 PetscBool dump = PETSC_FALSE, opsset, use_amat, missinginterpolate = PETSC_FALSE; 863 Mat dA, dB; 864 Vec tvec; 865 DM *dms; 866 PetscViewer viewer = NULL; 867 PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE; 868 PetscBool adaptInterpolation = mg->adaptInterpolation; 869 870 PetscFunctionBegin; 871 PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels with PCMGSetLevels() before setting up"); 872 n = mglevels[0]->levels; 873 /* FIX: Move this to PCSetFromOptions_MG? */ 874 if (mg->usedmfornumberoflevels) { 875 PetscInt levels; 876 PetscCall(DMGetRefineLevel(pc->dm, &levels)); 877 levels++; 878 if (levels > n) { /* the problem is now being solved on a finer grid */ 879 PetscCall(PCMGSetLevels(pc, levels, NULL)); 880 n = levels; 881 PetscCall(PCSetFromOptions(pc)); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 882 mglevels = mg->levels; 883 } 884 } 885 886 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 887 /* so use those from global PC */ 888 /* Is this what we always want? What if user wants to keep old one? */ 889 PetscCall(KSPGetOperatorsSet(mglevels[n - 1]->smoothd, NULL, &opsset)); 890 if (opsset) { 891 Mat mmat; 892 PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, NULL, &mmat)); 893 if (mmat == pc->pmat) opsset = PETSC_FALSE; 894 } 895 /* fine grid smoother inherits the reuse-pc flag */ 896 PetscCall(KSPGetPC(mglevels[n - 1]->smoothd, &cpc)); 897 cpc->reusepreconditioner = pc->reusepreconditioner; 898 PetscCall(KSPGetPC(mglevels[n - 1]->smoothu, &cpc)); 899 cpc->reusepreconditioner = pc->reusepreconditioner; 900 901 /* Create CR solvers */ 902 PetscCall(PCMGGetAdaptCR(pc, &doCR)); 903 if (doCR) { 904 const char *prefix; 905 906 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 907 for (i = 1; i < n; ++i) { 908 PC ipc, cr; 909 char crprefix[128]; 910 911 PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &mglevels[i]->cr)); 912 PetscCall(KSPSetNestLevel(mglevels[i]->cr, pc->kspnestlevel)); 913 PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE)); 914 PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->cr, (PetscObject)pc, n - i)); 915 PetscCall(KSPSetOptionsPrefix(mglevels[i]->cr, prefix)); 916 PetscCall(PetscObjectComposedDataSetInt((PetscObject)mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level)); 917 PetscCall(KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV)); 918 PetscCall(KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL)); 919 PetscCall(KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED)); 920 PetscCall(KSPGetPC(mglevels[i]->cr, &ipc)); 921 922 PetscCall(PCSetType(ipc, PCCOMPOSITE)); 923 PetscCall(PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE)); 924 PetscCall(PCCompositeAddPCType(ipc, PCSOR)); 925 PetscCall(CreateCR_Private(pc, i, &cr)); 926 PetscCall(PCCompositeAddPC(ipc, cr)); 927 PetscCall(PCDestroy(&cr)); 928 929 PetscCall(KSPSetTolerances(mglevels[i]->cr, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, mg->default_smoothd)); 930 PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 931 PetscCall(PetscSNPrintf(crprefix, 128, "mg_levels_%" PetscInt_FMT "_cr_", i)); 932 PetscCall(KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix)); 933 } 934 } 935 936 if (!opsset) { 937 PetscCall(PCGetUseAmat(pc, &use_amat)); 938 if (use_amat) { 939 PetscCall(PetscInfo(pc, "Using outer operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n")); 940 PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->mat, pc->pmat)); 941 } else { 942 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")); 943 PetscCall(KSPSetOperators(mglevels[n - 1]->smoothd, pc->pmat, pc->pmat)); 944 } 945 } 946 947 for (i = n - 1; i > 0; i--) { 948 if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 949 missinginterpolate = PETSC_TRUE; 950 break; 951 } 952 } 953 954 PetscCall(KSPGetOperators(mglevels[n - 1]->smoothd, &dA, &dB)); 955 if (dA == dB) dAeqdB = PETSC_TRUE; 956 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 958 if (pc->dm && !pc->setupcalled) { 959 /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 960 PetscCall(KSPSetDM(mglevels[n - 1]->smoothd, pc->dm)); 961 PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothd, PETSC_FALSE)); 962 if (mglevels[n - 1]->smoothd != mglevels[n - 1]->smoothu) { 963 PetscCall(KSPSetDM(mglevels[n - 1]->smoothu, pc->dm)); 964 PetscCall(KSPSetDMActive(mglevels[n - 1]->smoothu, PETSC_FALSE)); 965 } 966 if (mglevels[n - 1]->cr) { 967 PetscCall(KSPSetDM(mglevels[n - 1]->cr, pc->dm)); 968 PetscCall(KSPSetDMActive(mglevels[n - 1]->cr, PETSC_FALSE)); 969 } 970 } 971 972 /* 973 Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 974 Skipping for externally managed hierarchy (such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 975 */ 976 if (missinginterpolate && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 977 /* first see if we can compute a coarse space */ 978 if (mg->coarseSpaceType == PCMG_ADAPT_GDSW) { 979 for (i = n - 2; i > -1; i--) { 980 if (!mglevels[i + 1]->restrct && !mglevels[i + 1]->interpolate) { 981 PetscCall(PCMGComputeCoarseSpace_Internal(pc, i + 1, mg->coarseSpaceType, mg->Nc, NULL, &mglevels[i + 1]->coarseSpace)); 982 PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->coarseSpace)); 983 } 984 } 985 } else { /* construct the interpolation from the DMs */ 986 Mat p; 987 Vec rscale; 988 989 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); 990 PetscCall(PetscMalloc1(n, &dms)); 991 dms[n - 1] = pc->dm; 992 /* Separately create them so we do not get DMKSP interference between levels */ 993 for (i = n - 2; i > -1; i--) PetscCall(DMCoarsen(dms[i + 1], MPI_COMM_NULL, &dms[i])); 994 for (i = n - 2; i > -1; i--) { 995 DMKSP kdm; 996 PetscBool dmhasrestrict, dmhasinject; 997 998 PetscCall(KSPSetDM(mglevels[i]->smoothd, dms[i])); 999 if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothd, PETSC_FALSE)); 1000 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1001 PetscCall(KSPSetDM(mglevels[i]->smoothu, dms[i])); 1002 if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothu, PETSC_FALSE)); 1003 } 1004 if (mglevels[i]->cr) { 1005 PetscCall(KSPSetDM(mglevels[i]->cr, dms[i])); 1006 if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->cr, PETSC_FALSE)); 1007 } 1008 PetscCall(DMGetDMKSPWrite(dms[i], &kdm)); 1009 /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 1010 * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 1011 kdm->ops->computerhs = NULL; 1012 kdm->rhsctx = NULL; 1013 if (!mglevels[i + 1]->interpolate) { 1014 PetscCall(DMCreateInterpolation(dms[i], dms[i + 1], &p, &rscale)); 1015 PetscCall(PCMGSetInterpolation(pc, i + 1, p)); 1016 if (rscale) PetscCall(PCMGSetRScale(pc, i + 1, rscale)); 1017 PetscCall(VecDestroy(&rscale)); 1018 PetscCall(MatDestroy(&p)); 1019 } 1020 PetscCall(DMHasCreateRestriction(dms[i], &dmhasrestrict)); 1021 if (dmhasrestrict && !mglevels[i + 1]->restrct) { 1022 PetscCall(DMCreateRestriction(dms[i], dms[i + 1], &p)); 1023 PetscCall(PCMGSetRestriction(pc, i + 1, p)); 1024 PetscCall(MatDestroy(&p)); 1025 } 1026 PetscCall(DMHasCreateInjection(dms[i], &dmhasinject)); 1027 if (dmhasinject && !mglevels[i + 1]->inject) { 1028 PetscCall(DMCreateInjection(dms[i], dms[i + 1], &p)); 1029 PetscCall(PCMGSetInjection(pc, i + 1, p)); 1030 PetscCall(MatDestroy(&p)); 1031 } 1032 } 1033 1034 for (i = n - 2; i > -1; i--) PetscCall(DMDestroy(&dms[i])); 1035 PetscCall(PetscFree(dms)); 1036 } 1037 } 1038 1039 if (mg->galerkin < PC_MG_GALERKIN_NONE) { 1040 Mat A, B; 1041 PetscBool doA = PETSC_FALSE, doB = PETSC_FALSE; 1042 MatReuse reuse = MAT_INITIAL_MATRIX; 1043 1044 if (mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_BOTH) doB = PETSC_TRUE; 1045 if (mg->galerkin == PC_MG_GALERKIN_MAT || (mg->galerkin == PC_MG_GALERKIN_BOTH && dA != dB)) doA = PETSC_TRUE; 1046 if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 1047 for (i = n - 2; i > -1; i--) { 1048 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"); 1049 if (!mglevels[i + 1]->interpolate) PetscCall(PCMGSetInterpolation(pc, i + 1, mglevels[i + 1]->restrct)); 1050 if (!mglevels[i + 1]->restrct) PetscCall(PCMGSetRestriction(pc, i + 1, mglevels[i + 1]->interpolate)); 1051 if (reuse == MAT_REUSE_MATRIX) PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, &B)); 1052 if (doA) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dA, mglevels[i + 1]->interpolate, reuse, 1.0, &A)); 1053 if (doB) PetscCall(MatGalerkin(mglevels[i + 1]->restrct, dB, mglevels[i + 1]->interpolate, reuse, 1.0, &B)); 1054 /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 1055 if (!doA && dAeqdB) { 1056 if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)B)); 1057 A = B; 1058 } else if (!doA && reuse == MAT_INITIAL_MATRIX) { 1059 PetscCall(KSPGetOperators(mglevels[i]->smoothd, &A, NULL)); 1060 PetscCall(PetscObjectReference((PetscObject)A)); 1061 } 1062 if (!doB && dAeqdB) { 1063 if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)A)); 1064 B = A; 1065 } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 1066 PetscCall(KSPGetOperators(mglevels[i]->smoothd, NULL, &B)); 1067 PetscCall(PetscObjectReference((PetscObject)B)); 1068 } 1069 if (reuse == MAT_INITIAL_MATRIX) { 1070 PetscCall(KSPSetOperators(mglevels[i]->smoothd, A, B)); 1071 PetscCall(PetscObjectDereference((PetscObject)A)); 1072 PetscCall(PetscObjectDereference((PetscObject)B)); 1073 } 1074 dA = A; 1075 dB = B; 1076 } 1077 } 1078 1079 /* Adapt interpolation matrices */ 1080 if (adaptInterpolation) { 1081 for (i = 0; i < n; ++i) { 1082 if (!mglevels[i]->coarseSpace) PetscCall(PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i - 1]->coarseSpace, &mglevels[i]->coarseSpace)); 1083 if (i) PetscCall(PCMGAdaptInterpolator_Internal(pc, i, mglevels[i - 1]->smoothu, mglevels[i]->smoothu, mglevels[i - 1]->coarseSpace, mglevels[i]->coarseSpace)); 1084 } 1085 for (i = n - 2; i > -1; --i) PetscCall(PCMGRecomputeLevelOperators_Internal(pc, i)); 1086 } 1087 1088 if (needRestricts && pc->dm) { 1089 for (i = n - 2; i >= 0; i--) { 1090 DM dmfine, dmcoarse; 1091 Mat Restrict, Inject; 1092 Vec rscale; 1093 1094 PetscCall(KSPGetDM(mglevels[i + 1]->smoothd, &dmfine)); 1095 PetscCall(KSPGetDM(mglevels[i]->smoothd, &dmcoarse)); 1096 PetscCall(PCMGGetRestriction(pc, i + 1, &Restrict)); 1097 PetscCall(PCMGGetRScale(pc, i + 1, &rscale)); 1098 PetscCall(PCMGGetInjection(pc, i + 1, &Inject)); 1099 PetscCall(DMRestrict(dmfine, Restrict, rscale, Inject, dmcoarse)); 1100 } 1101 } 1102 1103 if (!pc->setupcalled) { 1104 for (i = 0; i < n; i++) PetscCall(KSPSetFromOptions(mglevels[i]->smoothd)); 1105 for (i = 1; i < n; i++) { 1106 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) PetscCall(KSPSetFromOptions(mglevels[i]->smoothu)); 1107 if (mglevels[i]->cr) PetscCall(KSPSetFromOptions(mglevels[i]->cr)); 1108 } 1109 /* insure that if either interpolation or restriction is set the other one is set */ 1110 for (i = 1; i < n; i++) { 1111 PetscCall(PCMGGetInterpolation(pc, i, NULL)); 1112 PetscCall(PCMGGetRestriction(pc, i, NULL)); 1113 } 1114 for (i = 0; i < n - 1; i++) { 1115 if (!mglevels[i]->b) { 1116 Vec *vec; 1117 PetscCall(KSPCreateVecs(mglevels[i]->smoothd, 1, &vec, 0, NULL)); 1118 PetscCall(PCMGSetRhs(pc, i, *vec)); 1119 PetscCall(VecDestroy(vec)); 1120 PetscCall(PetscFree(vec)); 1121 } 1122 if (!mglevels[i]->r && i) { 1123 PetscCall(VecDuplicate(mglevels[i]->b, &tvec)); 1124 PetscCall(PCMGSetR(pc, i, tvec)); 1125 PetscCall(VecDestroy(&tvec)); 1126 } 1127 if (!mglevels[i]->x) { 1128 PetscCall(VecDuplicate(mglevels[i]->b, &tvec)); 1129 PetscCall(PCMGSetX(pc, i, tvec)); 1130 PetscCall(VecDestroy(&tvec)); 1131 } 1132 if (doCR) { 1133 PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crx)); 1134 PetscCall(VecDuplicate(mglevels[i]->b, &mglevels[i]->crb)); 1135 PetscCall(VecZeroEntries(mglevels[i]->crb)); 1136 } 1137 } 1138 if (n != 1 && !mglevels[n - 1]->r) { 1139 /* PCMGSetR() on the finest level if user did not supply it */ 1140 Vec *vec; 1141 1142 PetscCall(KSPCreateVecs(mglevels[n - 1]->smoothd, 1, &vec, 0, NULL)); 1143 PetscCall(PCMGSetR(pc, n - 1, *vec)); 1144 PetscCall(VecDestroy(vec)); 1145 PetscCall(PetscFree(vec)); 1146 } 1147 if (doCR) { 1148 PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crx)); 1149 PetscCall(VecDuplicate(mglevels[n - 1]->r, &mglevels[n - 1]->crb)); 1150 PetscCall(VecZeroEntries(mglevels[n - 1]->crb)); 1151 } 1152 } 1153 1154 if (pc->dm) { 1155 /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 1156 for (i = 0; i < n - 1; i++) { 1157 if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 1158 } 1159 } 1160 // We got here (PCSetUp_MG) because the matrix has changed, which means the smoother needs to be set up again (e.g., 1161 // new diagonal for Jacobi). Setting it here allows it to be logged under PCSetUp rather than deep inside a PCApply. 1162 if (mglevels[n - 1]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[n - 1]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 1163 1164 for (i = 1; i < n; i++) { 1165 if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1) { 1166 /* if doing only down then initial guess is zero */ 1167 PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd, PETSC_TRUE)); 1168 } 1169 if (mglevels[i]->cr) PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 1170 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 1171 PetscCall(KSPSetUp(mglevels[i]->smoothd)); 1172 if (mglevels[i]->smoothd->reason) pc->failedreason = PC_SUBPC_ERROR; 1173 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 1174 if (!mglevels[i]->residual) { 1175 Mat mat; 1176 1177 PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL)); 1178 PetscCall(PCMGSetResidual(pc, i, PCMGResidualDefault, mat)); 1179 } 1180 if (!mglevels[i]->residualtranspose) { 1181 Mat mat; 1182 1183 PetscCall(KSPGetOperators(mglevels[i]->smoothd, &mat, NULL)); 1184 PetscCall(PCMGSetResidualTranspose(pc, i, PCMGResidualTransposeDefault, mat)); 1185 } 1186 } 1187 for (i = 1; i < n; i++) { 1188 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 1189 Mat downmat, downpmat; 1190 1191 /* check if operators have been set for up, if not use down operators to set them */ 1192 PetscCall(KSPGetOperatorsSet(mglevels[i]->smoothu, &opsset, NULL)); 1193 if (!opsset) { 1194 PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat)); 1195 PetscCall(KSPSetOperators(mglevels[i]->smoothu, downmat, downpmat)); 1196 } 1197 1198 PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothu, PETSC_TRUE)); 1199 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 1200 PetscCall(KSPSetUp(mglevels[i]->smoothu)); 1201 if (mglevels[i]->smoothu->reason) pc->failedreason = PC_SUBPC_ERROR; 1202 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 1203 } 1204 if (mglevels[i]->cr) { 1205 Mat downmat, downpmat; 1206 1207 /* check if operators have been set for up, if not use down operators to set them */ 1208 PetscCall(KSPGetOperatorsSet(mglevels[i]->cr, &opsset, NULL)); 1209 if (!opsset) { 1210 PetscCall(KSPGetOperators(mglevels[i]->smoothd, &downmat, &downpmat)); 1211 PetscCall(KSPSetOperators(mglevels[i]->cr, downmat, downpmat)); 1212 } 1213 1214 PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 1215 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 1216 PetscCall(KSPSetUp(mglevels[i]->cr)); 1217 if (mglevels[i]->cr->reason) pc->failedreason = PC_SUBPC_ERROR; 1218 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup, 0, 0, 0, 0)); 1219 } 1220 } 1221 1222 if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0)); 1223 PetscCall(KSPSetUp(mglevels[0]->smoothd)); 1224 if (mglevels[0]->smoothd->reason) pc->failedreason = PC_SUBPC_ERROR; 1225 if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[0]->eventsmoothsetup, 0, 0, 0, 0)); 1226 1227 /* 1228 Dump the interpolation/restriction matrices plus the 1229 Jacobian/stiffness on each level. This allows MATLAB users to 1230 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 1231 1232 Only support one or the other at the same time. 1233 */ 1234 #if defined(PETSC_USE_SOCKET_VIEWER) 1235 PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_matlab", &dump, NULL)); 1236 if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 1237 dump = PETSC_FALSE; 1238 #endif 1239 PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options, ((PetscObject)pc)->prefix, "-pc_mg_dump_binary", &dump, NULL)); 1240 if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 1241 1242 if (viewer) { 1243 for (i = 1; i < n; i++) PetscCall(MatView(mglevels[i]->restrct, viewer)); 1244 for (i = 0; i < n; i++) { 1245 PetscCall(KSPGetPC(mglevels[i]->smoothd, &pc)); 1246 PetscCall(MatView(pc->mat, viewer)); 1247 } 1248 } 1249 PetscFunctionReturn(PETSC_SUCCESS); 1250 } 1251 1252 PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) 1253 { 1254 PC_MG *mg = (PC_MG *)pc->data; 1255 1256 PetscFunctionBegin; 1257 *levels = mg->nlevels; 1258 PetscFunctionReturn(PETSC_SUCCESS); 1259 } 1260 1261 /*@ 1262 PCMGGetLevels - Gets the number of levels to use with `PCMG`. 1263 1264 Not Collective 1265 1266 Input Parameter: 1267 . pc - the preconditioner context 1268 1269 Output Parameter: 1270 . levels - the number of levels 1271 1272 Level: advanced 1273 1274 .seealso: [](ch_ksp), `PCMG`, `PCMGSetLevels()` 1275 @*/ 1276 PetscErrorCode PCMGGetLevels(PC pc, PetscInt *levels) 1277 { 1278 PetscFunctionBegin; 1279 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1280 PetscAssertPointer(levels, 2); 1281 *levels = 0; 1282 PetscTryMethod(pc, "PCMGGetLevels_C", (PC, PetscInt *), (pc, levels)); 1283 PetscFunctionReturn(PETSC_SUCCESS); 1284 } 1285 1286 /*@ 1287 PCMGGetGridComplexity - compute operator and grid complexity of the `PCMG` hierarchy 1288 1289 Input Parameter: 1290 . pc - the preconditioner context 1291 1292 Output Parameters: 1293 + gc - grid complexity = sum_i(n_i) / n_0 1294 - oc - operator complexity = sum_i(nnz_i) / nnz_0 1295 1296 Level: advanced 1297 1298 Note: 1299 This is often call the operator complexity in multigrid literature 1300 1301 .seealso: [](ch_ksp), `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()` 1302 @*/ 1303 PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc, PetscReal *oc) 1304 { 1305 PC_MG *mg = (PC_MG *)pc->data; 1306 PC_MG_Levels **mglevels = mg->levels; 1307 PetscInt lev, N; 1308 PetscLogDouble nnz0 = 0, sgc = 0, soc = 0, n0 = 0; 1309 MatInfo info; 1310 1311 PetscFunctionBegin; 1312 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1313 if (gc) PetscAssertPointer(gc, 2); 1314 if (oc) PetscAssertPointer(oc, 3); 1315 if (!pc->setupcalled) { 1316 if (gc) *gc = 0; 1317 if (oc) *oc = 0; 1318 PetscFunctionReturn(PETSC_SUCCESS); 1319 } 1320 PetscCheck(mg->nlevels > 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "MG has no levels"); 1321 for (lev = 0; lev < mg->nlevels; lev++) { 1322 Mat dB; 1323 PetscCall(KSPGetOperators(mglevels[lev]->smoothd, NULL, &dB)); 1324 PetscCall(MatGetInfo(dB, MAT_GLOBAL_SUM, &info)); /* global reduction */ 1325 PetscCall(MatGetSize(dB, &N, NULL)); 1326 sgc += N; 1327 soc += info.nz_used; 1328 if (lev == mg->nlevels - 1) { 1329 nnz0 = info.nz_used; 1330 n0 = N; 1331 } 1332 } 1333 PetscCheck(n0 > 0 && gc, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number for grid points on finest level is not available"); 1334 *gc = (PetscReal)(sgc / n0); 1335 if (nnz0 > 0 && oc) *oc = (PetscReal)(soc / nnz0); 1336 PetscFunctionReturn(PETSC_SUCCESS); 1337 } 1338 1339 /*@ 1340 PCMGSetType - Determines the form of multigrid to use, either 1341 multiplicative, additive, full, or the Kaskade algorithm. 1342 1343 Logically Collective 1344 1345 Input Parameters: 1346 + pc - the preconditioner context 1347 - form - multigrid form, one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE` 1348 1349 Options Database Key: 1350 . -pc_mg_type <form> - Sets <form>, one of multiplicative, additive, full, kaskade 1351 1352 Level: advanced 1353 1354 .seealso: [](ch_ksp), `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGGetType()`, `PCMGCycleType` 1355 @*/ 1356 PetscErrorCode PCMGSetType(PC pc, PCMGType form) 1357 { 1358 PC_MG *mg = (PC_MG *)pc->data; 1359 1360 PetscFunctionBegin; 1361 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1362 PetscValidLogicalCollectiveEnum(pc, form, 2); 1363 mg->am = form; 1364 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 1365 else pc->ops->applyrichardson = NULL; 1366 PetscFunctionReturn(PETSC_SUCCESS); 1367 } 1368 1369 /*@ 1370 PCMGGetType - Finds the form of multigrid the `PCMG` is using multiplicative, additive, full, or the Kaskade algorithm. 1371 1372 Logically Collective 1373 1374 Input Parameter: 1375 . pc - the preconditioner context 1376 1377 Output Parameter: 1378 . type - one of `PC_MG_MULTIPLICATIVE`, `PC_MG_ADDITIVE`, `PC_MG_FULL`, `PC_MG_KASKADE`, `PCMGCycleType` 1379 1380 Level: advanced 1381 1382 .seealso: [](ch_ksp), `PCMGType`, `PCMG`, `PCMGGetLevels()`, `PCMGSetLevels()`, `PCMGSetType()` 1383 @*/ 1384 PetscErrorCode PCMGGetType(PC pc, PCMGType *type) 1385 { 1386 PC_MG *mg = (PC_MG *)pc->data; 1387 1388 PetscFunctionBegin; 1389 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1390 *type = mg->am; 1391 PetscFunctionReturn(PETSC_SUCCESS); 1392 } 1393 1394 /*@ 1395 PCMGSetCycleType - Sets the type cycles to use. Use `PCMGSetCycleTypeOnLevel()` for more 1396 complicated cycling. 1397 1398 Logically Collective 1399 1400 Input Parameters: 1401 + pc - the multigrid context 1402 - n - either `PC_MG_CYCLE_V` or `PC_MG_CYCLE_W` 1403 1404 Options Database Key: 1405 . -pc_mg_cycle_type <v,w> - provide the cycle desired 1406 1407 Level: advanced 1408 1409 .seealso: [](ch_ksp), `PCMG`, `PCMGSetCycleTypeOnLevel()`, `PCMGType`, `PCMGCycleType` 1410 @*/ 1411 PetscErrorCode PCMGSetCycleType(PC pc, PCMGCycleType n) 1412 { 1413 PC_MG *mg = (PC_MG *)pc->data; 1414 PC_MG_Levels **mglevels = mg->levels; 1415 PetscInt i, levels; 1416 1417 PetscFunctionBegin; 1418 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1419 PetscValidLogicalCollectiveEnum(pc, n, 2); 1420 PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling"); 1421 levels = mglevels[0]->levels; 1422 for (i = 0; i < levels; i++) mglevels[i]->cycles = n; 1423 PetscFunctionReturn(PETSC_SUCCESS); 1424 } 1425 1426 /*@ 1427 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 1428 of multigrid when `PCMGType` is `PC_MG_MULTIPLICATIVE` 1429 1430 Logically Collective 1431 1432 Input Parameters: 1433 + pc - the multigrid context 1434 - n - number of cycles (default is 1) 1435 1436 Options Database Key: 1437 . -pc_mg_multiplicative_cycles n - set the number of cycles 1438 1439 Level: advanced 1440 1441 Note: 1442 This is not associated with setting a v or w cycle, that is set with `PCMGSetCycleType()` 1443 1444 .seealso: [](ch_ksp), `PCMGSetCycleTypeOnLevel()`, `PCMGSetCycleType()`, `PCMGCycleType`, `PCMGType` 1445 @*/ 1446 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc, PetscInt n) 1447 { 1448 PC_MG *mg = (PC_MG *)pc->data; 1449 1450 PetscFunctionBegin; 1451 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1452 PetscValidLogicalCollectiveInt(pc, n, 2); 1453 mg->cyclesperpcapply = n; 1454 PetscFunctionReturn(PETSC_SUCCESS); 1455 } 1456 1457 /* 1458 Since the finest level KSP shares the original matrix (of the entire system), it's preconditioner 1459 should not be updated if the whole PC is supposed to reuse the preconditioner 1460 */ 1461 static PetscErrorCode PCSetReusePreconditioner_MG(PC pc, PetscBool flag) 1462 { 1463 PC_MG *mg = (PC_MG *)pc->data; 1464 PC_MG_Levels **mglevels = mg->levels; 1465 PetscInt levels; 1466 PC tpc; 1467 1468 PetscFunctionBegin; 1469 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1470 PetscValidLogicalCollectiveBool(pc, flag, 2); 1471 if (mglevels) { 1472 levels = mglevels[0]->levels; 1473 PetscCall(KSPGetPC(mglevels[levels - 1]->smoothd, &tpc)); 1474 tpc->reusepreconditioner = flag; 1475 PetscCall(KSPGetPC(mglevels[levels - 1]->smoothu, &tpc)); 1476 tpc->reusepreconditioner = flag; 1477 } 1478 PetscFunctionReturn(PETSC_SUCCESS); 1479 } 1480 1481 static PetscErrorCode PCMGSetGalerkin_MG(PC pc, PCMGGalerkinType use) 1482 { 1483 PC_MG *mg = (PC_MG *)pc->data; 1484 1485 PetscFunctionBegin; 1486 mg->galerkin = use; 1487 PetscFunctionReturn(PETSC_SUCCESS); 1488 } 1489 1490 /*@ 1491 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 1492 finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1493 1494 Logically Collective 1495 1496 Input Parameters: 1497 + pc - the multigrid context 1498 - use - one of `PC_MG_GALERKIN_BOTH`, `PC_MG_GALERKIN_PMAT`, `PC_MG_GALERKIN_MAT`, or `PC_MG_GALERKIN_NONE` 1499 1500 Options Database Key: 1501 . -pc_mg_galerkin <both,pmat,mat,none> - set the matrices to form via the Galerkin process 1502 1503 Level: intermediate 1504 1505 Note: 1506 Some codes that use `PCMG` such as `PCGAMG` use Galerkin internally while constructing the hierarchy and thus do not 1507 use the `PCMG` construction of the coarser grids. 1508 1509 .seealso: [](ch_ksp), `PCMG`, `PCMGGetGalerkin()`, `PCMGGalerkinType` 1510 @*/ 1511 PetscErrorCode PCMGSetGalerkin(PC pc, PCMGGalerkinType use) 1512 { 1513 PetscFunctionBegin; 1514 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1515 PetscTryMethod(pc, "PCMGSetGalerkin_C", (PC, PCMGGalerkinType), (pc, use)); 1516 PetscFunctionReturn(PETSC_SUCCESS); 1517 } 1518 1519 /*@ 1520 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. A_i-1 = r_i * A_i * p_i 1521 1522 Not Collective 1523 1524 Input Parameter: 1525 . pc - the multigrid context 1526 1527 Output Parameter: 1528 . galerkin - one of `PC_MG_GALERKIN_BOTH`,`PC_MG_GALERKIN_PMAT`,`PC_MG_GALERKIN_MAT`, `PC_MG_GALERKIN_NONE`, or `PC_MG_GALERKIN_EXTERNAL` 1529 1530 Level: intermediate 1531 1532 .seealso: [](ch_ksp), `PCMG`, `PCMGSetGalerkin()`, `PCMGGalerkinType` 1533 @*/ 1534 PetscErrorCode PCMGGetGalerkin(PC pc, PCMGGalerkinType *galerkin) 1535 { 1536 PC_MG *mg = (PC_MG *)pc->data; 1537 1538 PetscFunctionBegin; 1539 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1540 *galerkin = mg->galerkin; 1541 PetscFunctionReturn(PETSC_SUCCESS); 1542 } 1543 1544 static PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt) 1545 { 1546 PC_MG *mg = (PC_MG *)pc->data; 1547 1548 PetscFunctionBegin; 1549 mg->adaptInterpolation = adapt; 1550 PetscFunctionReturn(PETSC_SUCCESS); 1551 } 1552 1553 static PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt) 1554 { 1555 PC_MG *mg = (PC_MG *)pc->data; 1556 1557 PetscFunctionBegin; 1558 *adapt = mg->adaptInterpolation; 1559 PetscFunctionReturn(PETSC_SUCCESS); 1560 } 1561 1562 static PetscErrorCode PCMGSetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType ctype) 1563 { 1564 PC_MG *mg = (PC_MG *)pc->data; 1565 1566 PetscFunctionBegin; 1567 mg->adaptInterpolation = ctype != PCMG_ADAPT_NONE ? PETSC_TRUE : PETSC_FALSE; 1568 mg->coarseSpaceType = ctype; 1569 PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_BOTH)); 1570 PetscFunctionReturn(PETSC_SUCCESS); 1571 } 1572 1573 static PetscErrorCode PCMGGetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType *ctype) 1574 { 1575 PC_MG *mg = (PC_MG *)pc->data; 1576 1577 PetscFunctionBegin; 1578 *ctype = mg->coarseSpaceType; 1579 PetscFunctionReturn(PETSC_SUCCESS); 1580 } 1581 1582 static PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr) 1583 { 1584 PC_MG *mg = (PC_MG *)pc->data; 1585 1586 PetscFunctionBegin; 1587 mg->compatibleRelaxation = cr; 1588 PetscFunctionReturn(PETSC_SUCCESS); 1589 } 1590 1591 static PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr) 1592 { 1593 PC_MG *mg = (PC_MG *)pc->data; 1594 1595 PetscFunctionBegin; 1596 *cr = mg->compatibleRelaxation; 1597 PetscFunctionReturn(PETSC_SUCCESS); 1598 } 1599 1600 /*@ 1601 PCMGSetAdaptCoarseSpaceType - Set the type of adaptive coarse space. 1602 1603 Adapts or creates the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated. 1604 1605 Logically Collective 1606 1607 Input Parameters: 1608 + pc - the multigrid context 1609 - ctype - the type of coarse space 1610 1611 Options Database Keys: 1612 + -pc_mg_adapt_interp_n <int> - The number of modes to use 1613 - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: none, `polynomial`, `harmonic`, `eigenvector`, `generalized_eigenvector`, `gdsw` 1614 1615 Level: intermediate 1616 1617 Note: 1618 Requires a `DM` with specific functionality be attached to the `PC`. 1619 1620 .seealso: [](ch_ksp), `PCMG`, `PCMGCoarseSpaceType`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()`, `DM` 1621 @*/ 1622 PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType ctype) 1623 { 1624 PetscFunctionBegin; 1625 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1626 PetscValidLogicalCollectiveEnum(pc, ctype, 2); 1627 PetscTryMethod(pc, "PCMGSetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType), (pc, ctype)); 1628 PetscFunctionReturn(PETSC_SUCCESS); 1629 } 1630 1631 /*@ 1632 PCMGGetAdaptCoarseSpaceType - Get the type of adaptive coarse space. 1633 1634 Not Collective 1635 1636 Input Parameter: 1637 . pc - the multigrid context 1638 1639 Output Parameter: 1640 . ctype - the type of coarse space 1641 1642 Level: intermediate 1643 1644 .seealso: [](ch_ksp), `PCMG`, `PCMGCoarseSpaceType`, `PCMGSetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()`, `PCMGSetAdaptInterpolation()` 1645 @*/ 1646 PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType *ctype) 1647 { 1648 PetscFunctionBegin; 1649 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1650 PetscAssertPointer(ctype, 2); 1651 PetscUseMethod(pc, "PCMGGetAdaptCoarseSpaceType_C", (PC, PCMGCoarseSpaceType *), (pc, ctype)); 1652 PetscFunctionReturn(PETSC_SUCCESS); 1653 } 1654 1655 /* MATT: REMOVE? */ 1656 /*@ 1657 PCMGSetAdaptInterpolation - Adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated. 1658 1659 Logically Collective 1660 1661 Input Parameters: 1662 + pc - the multigrid context 1663 - adapt - flag for adaptation of the interpolator 1664 1665 Options Database Keys: 1666 + -pc_mg_adapt_interp - Turn on adaptation 1667 . -pc_mg_adapt_interp_n <int> - The number of modes to use, should be divisible by dimension 1668 - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector 1669 1670 Level: intermediate 1671 1672 .seealso: [](ch_ksp), `PCMG`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 1673 @*/ 1674 PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt) 1675 { 1676 PetscFunctionBegin; 1677 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1678 PetscTryMethod(pc, "PCMGSetAdaptInterpolation_C", (PC, PetscBool), (pc, adapt)); 1679 PetscFunctionReturn(PETSC_SUCCESS); 1680 } 1681 1682 /*@ 1683 PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, 1684 and thus accurately interpolated. 1685 1686 Not Collective 1687 1688 Input Parameter: 1689 . pc - the multigrid context 1690 1691 Output Parameter: 1692 . adapt - flag for adaptation of the interpolator 1693 1694 Level: intermediate 1695 1696 .seealso: [](ch_ksp), `PCMG`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 1697 @*/ 1698 PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt) 1699 { 1700 PetscFunctionBegin; 1701 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1702 PetscAssertPointer(adapt, 2); 1703 PetscUseMethod(pc, "PCMGGetAdaptInterpolation_C", (PC, PetscBool *), (pc, adapt)); 1704 PetscFunctionReturn(PETSC_SUCCESS); 1705 } 1706 1707 /*@ 1708 PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation. 1709 1710 Logically Collective 1711 1712 Input Parameters: 1713 + pc - the multigrid context 1714 - cr - flag for compatible relaxation 1715 1716 Options Database Key: 1717 . -pc_mg_adapt_cr - Turn on compatible relaxation 1718 1719 Level: intermediate 1720 1721 .seealso: [](ch_ksp), `PCMG`, `PCMGGetAdaptCR()`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 1722 @*/ 1723 PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr) 1724 { 1725 PetscFunctionBegin; 1726 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1727 PetscTryMethod(pc, "PCMGSetAdaptCR_C", (PC, PetscBool), (pc, cr)); 1728 PetscFunctionReturn(PETSC_SUCCESS); 1729 } 1730 1731 /*@ 1732 PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation. 1733 1734 Not Collective 1735 1736 Input Parameter: 1737 . pc - the multigrid context 1738 1739 Output Parameter: 1740 . cr - flag for compatible relaxaion 1741 1742 Level: intermediate 1743 1744 .seealso: [](ch_ksp), `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 1745 @*/ 1746 PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr) 1747 { 1748 PetscFunctionBegin; 1749 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1750 PetscAssertPointer(cr, 2); 1751 PetscUseMethod(pc, "PCMGGetAdaptCR_C", (PC, PetscBool *), (pc, cr)); 1752 PetscFunctionReturn(PETSC_SUCCESS); 1753 } 1754 1755 /*@ 1756 PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1757 on all levels. Use `PCMGDistinctSmoothUp()` to create separate up and down smoothers if you want different numbers of 1758 pre- and post-smoothing steps. 1759 1760 Logically Collective 1761 1762 Input Parameters: 1763 + pc - the multigrid context 1764 - n - the number of smoothing steps 1765 1766 Options Database Key: 1767 . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 1768 1769 Level: advanced 1770 1771 Note: 1772 This does not set a value on the coarsest grid, since we assume that there is no separate smooth up on the coarsest grid. 1773 1774 .seealso: [](ch_ksp), `PCMG`, `PCMGSetDistinctSmoothUp()` 1775 @*/ 1776 PetscErrorCode PCMGSetNumberSmooth(PC pc, PetscInt n) 1777 { 1778 PC_MG *mg = (PC_MG *)pc->data; 1779 PC_MG_Levels **mglevels = mg->levels; 1780 PetscInt i, levels; 1781 1782 PetscFunctionBegin; 1783 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1784 PetscValidLogicalCollectiveInt(pc, n, 2); 1785 PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling"); 1786 levels = mglevels[0]->levels; 1787 1788 for (i = 1; i < levels; i++) { 1789 PetscCall(KSPSetTolerances(mglevels[i]->smoothu, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, n)); 1790 PetscCall(KSPSetTolerances(mglevels[i]->smoothd, PETSC_CURRENT, PETSC_CURRENT, PETSC_CURRENT, n)); 1791 mg->default_smoothu = n; 1792 mg->default_smoothd = n; 1793 } 1794 PetscFunctionReturn(PETSC_SUCCESS); 1795 } 1796 1797 /*@ 1798 PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate `KSP` from the down (pre) smoother on all levels 1799 and adds the suffix _up to the options name 1800 1801 Logically Collective 1802 1803 Input Parameter: 1804 . pc - the preconditioner context 1805 1806 Options Database Key: 1807 . -pc_mg_distinct_smoothup <bool> - use distinct smoothing objects 1808 1809 Level: advanced 1810 1811 Note: 1812 This does not set a value on the coarsest grid, since we assume that there is no separate smooth up on the coarsest grid. 1813 1814 .seealso: [](ch_ksp), `PCMG`, `PCMGSetNumberSmooth()` 1815 @*/ 1816 PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1817 { 1818 PC_MG *mg = (PC_MG *)pc->data; 1819 PC_MG_Levels **mglevels = mg->levels; 1820 PetscInt i, levels; 1821 KSP subksp; 1822 1823 PetscFunctionBegin; 1824 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1825 PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ORDER, "Must set MG levels with PCMGSetLevels() before calling"); 1826 levels = mglevels[0]->levels; 1827 1828 for (i = 1; i < levels; i++) { 1829 const char *prefix = NULL; 1830 /* make sure smoother up and down are different */ 1831 PetscCall(PCMGGetSmootherUp(pc, i, &subksp)); 1832 PetscCall(KSPGetOptionsPrefix(mglevels[i]->smoothd, &prefix)); 1833 PetscCall(KSPSetOptionsPrefix(subksp, prefix)); 1834 PetscCall(KSPAppendOptionsPrefix(subksp, "up_")); 1835 } 1836 PetscFunctionReturn(PETSC_SUCCESS); 1837 } 1838 1839 /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 1840 static PetscErrorCode PCGetInterpolations_MG(PC pc, PetscInt *num_levels, Mat *interpolations[]) 1841 { 1842 PC_MG *mg = (PC_MG *)pc->data; 1843 PC_MG_Levels **mglevels = mg->levels; 1844 Mat *mat; 1845 PetscInt l; 1846 1847 PetscFunctionBegin; 1848 PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 1849 PetscCall(PetscMalloc1(mg->nlevels, &mat)); 1850 for (l = 1; l < mg->nlevels; l++) { 1851 mat[l - 1] = mglevels[l]->interpolate; 1852 PetscCall(PetscObjectReference((PetscObject)mat[l - 1])); 1853 } 1854 *num_levels = mg->nlevels; 1855 *interpolations = mat; 1856 PetscFunctionReturn(PETSC_SUCCESS); 1857 } 1858 1859 /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 1860 static PetscErrorCode PCGetCoarseOperators_MG(PC pc, PetscInt *num_levels, Mat *coarseOperators[]) 1861 { 1862 PC_MG *mg = (PC_MG *)pc->data; 1863 PC_MG_Levels **mglevels = mg->levels; 1864 PetscInt l; 1865 Mat *mat; 1866 1867 PetscFunctionBegin; 1868 PetscCheck(mglevels, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set MG levels before calling"); 1869 PetscCall(PetscMalloc1(mg->nlevels, &mat)); 1870 for (l = 0; l < mg->nlevels - 1; l++) { 1871 PetscCall(KSPGetOperators(mglevels[l]->smoothd, NULL, &mat[l])); 1872 PetscCall(PetscObjectReference((PetscObject)mat[l])); 1873 } 1874 *num_levels = mg->nlevels; 1875 *coarseOperators = mat; 1876 PetscFunctionReturn(PETSC_SUCCESS); 1877 } 1878 1879 /*@C 1880 PCMGRegisterCoarseSpaceConstructor - Adds a method to the `PCMG` package for coarse space construction. 1881 1882 Not Collective, No Fortran Support 1883 1884 Input Parameters: 1885 + name - name of the constructor 1886 - function - constructor routine, see `PCMGCoarseSpaceConstructorFn` 1887 1888 Level: advanced 1889 1890 Developer Notes: 1891 This does not appear to be used anywhere 1892 1893 .seealso: [](ch_ksp), `PCMGCoarseSpaceConstructorFn`, `PCMG`, `PCMGGetCoarseSpaceConstructor()`, `PCRegister()` 1894 @*/ 1895 PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PCMGCoarseSpaceConstructorFn *function) 1896 { 1897 PetscFunctionBegin; 1898 PetscCall(PCInitializePackage()); 1899 PetscCall(PetscFunctionListAdd(&PCMGCoarseList, name, function)); 1900 PetscFunctionReturn(PETSC_SUCCESS); 1901 } 1902 1903 /*@C 1904 PCMGGetCoarseSpaceConstructor - Returns the given coarse space construction method. 1905 1906 Not Collective, No Fortran Support 1907 1908 Input Parameter: 1909 . name - name of the constructor 1910 1911 Output Parameter: 1912 . function - constructor routine 1913 1914 Level: advanced 1915 1916 .seealso: [](ch_ksp), `PCMGCoarseSpaceConstructorFn`, `PCMG`, `PCMGRegisterCoarseSpaceConstructor()`, `PCRegister()` 1917 @*/ 1918 PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PCMGCoarseSpaceConstructorFn **function) 1919 { 1920 PetscFunctionBegin; 1921 PetscCall(PetscFunctionListFind(PCMGCoarseList, name, function)); 1922 PetscFunctionReturn(PETSC_SUCCESS); 1923 } 1924 1925 /*MC 1926 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 1927 information about the coarser grid matrices and restriction/interpolation operators. 1928 1929 Options Database Keys: 1930 + -pc_mg_levels <nlevels> - number of levels including finest 1931 . -pc_mg_cycle_type <v,w> - provide the cycle desired 1932 . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 1933 . -pc_mg_log - log information about time spent on each level of the solver 1934 . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 1935 . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 1936 . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 1937 . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1938 to a `PETSCVIEWERSOCKET` for reading from MATLAB. 1939 - -pc_mg_dump_binary -dumps the matrices for each level and the restriction/interpolation matrices 1940 to the binary output file called binaryoutput 1941 1942 Level: intermediate 1943 1944 Notes: 1945 The Krylov solver (if any) and preconditioner (smoother) and their parameters are controlled from the options database with the standard 1946 options database keywords prefixed with `-mg_levels_` to affect all the levels but the coarsest, which is controlled with `-mg_coarse_`, 1947 and the finest where `-mg_fine_` can override `-mg_levels_`. One can set different preconditioners etc on specific levels with the prefix 1948 `-mg_levels_n_` where `n` is the level number (zero being the coarse level. For example 1949 .vb 1950 -mg_levels_ksp_type gmres -mg_levels_pc_type bjacobi -mg_coarse_pc_type svd -mg_levels_2_pc_type sor 1951 .ve 1952 These options also work for controlling the smoothers etc inside `PCGAMG` 1953 1954 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 1955 1956 When run with a single level the smoother options are used on that level NOT the coarse grid solver options 1957 1958 When run with `KSPRICHARDSON` the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 1959 is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 1960 (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 1961 residual is computed at the end of each cycle. 1962 1963 .seealso: [](sec_mg), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCEXOTIC`, `PCGAMG`, `PCML`, `PCHYPRE` 1964 `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `PCMGSetCycleType()`, 1965 `PCMGSetDistinctSmoothUp()`, `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`, 1966 `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, 1967 `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()`, 1968 `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetAdaptCoarseSpaceType()` 1969 M*/ 1970 1971 PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 1972 { 1973 PC_MG *mg; 1974 1975 PetscFunctionBegin; 1976 PetscCall(PetscNew(&mg)); 1977 pc->data = mg; 1978 mg->nlevels = -1; 1979 mg->am = PC_MG_MULTIPLICATIVE; 1980 mg->galerkin = PC_MG_GALERKIN_NONE; 1981 mg->adaptInterpolation = PETSC_FALSE; 1982 mg->Nc = -1; 1983 mg->eigenvalue = -1; 1984 1985 pc->useAmat = PETSC_TRUE; 1986 1987 pc->ops->apply = PCApply_MG; 1988 pc->ops->applytranspose = PCApplyTranspose_MG; 1989 pc->ops->matapply = PCMatApply_MG; 1990 pc->ops->matapplytranspose = PCMatApplyTranspose_MG; 1991 pc->ops->setup = PCSetUp_MG; 1992 pc->ops->reset = PCReset_MG; 1993 pc->ops->destroy = PCDestroy_MG; 1994 pc->ops->setfromoptions = PCSetFromOptions_MG; 1995 pc->ops->view = PCView_MG; 1996 1997 PetscCall(PetscObjectComposedDataRegister(&mg->eigenvalue)); 1998 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetGalerkin_C", PCMGSetGalerkin_MG)); 1999 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSetReusePreconditioner_C", PCSetReusePreconditioner_MG)); 2000 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetLevels_C", PCMGGetLevels_MG)); 2001 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetLevels_C", PCMGSetLevels_MG)); 2002 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetInterpolations_C", PCGetInterpolations_MG)); 2003 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGetCoarseOperators_C", PCGetCoarseOperators_MG)); 2004 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptInterpolation_C", PCMGSetAdaptInterpolation_MG)); 2005 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptInterpolation_C", PCMGGetAdaptInterpolation_MG)); 2006 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCR_C", PCMGSetAdaptCR_MG)); 2007 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCR_C", PCMGGetAdaptCR_MG)); 2008 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGSetAdaptCoarseSpaceType_C", PCMGSetAdaptCoarseSpaceType_MG)); 2009 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMGGetAdaptCoarseSpaceType_C", PCMGGetAdaptCoarseSpaceType_MG)); 2010 PetscFunctionReturn(PETSC_SUCCESS); 2011 } 2012