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