1 #include <petscdm.h> 2 #include <petsc/private/hashmapi.h> 3 #include <petsc/private/matimpl.h> 4 #include <petsc/private/pcmgimpl.h> 5 #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 6 7 typedef struct { 8 PC innerpc; /* A MG inner PC (Hypre or PCGAMG) to setup interpolations and coarse operators */ 9 char *innerpctype; /* PCGAMG or PCHYPRE */ 10 PetscBool reuseinterp; /* A flag indicates if or not to reuse the interpolations */ 11 PetscBool subcoarsening; /* If or not to use a subspace-based coarsening algorithm */ 12 PetscBool usematmaij; /* If or not to use MatMAIJ for saving memory */ 13 PetscInt component; /* Which subspace is used for the subspace-based coarsening algorithm? */ 14 } PC_HMG; 15 16 PetscErrorCode PCReset_MG(PC); 17 18 static PetscErrorCode PCHMGExtractSubMatrix_Private(Mat pmat, Mat *submat, MatReuse reuse, PetscInt component, PetscInt blocksize) 19 { 20 IS isrow; 21 PetscInt rstart, rend; 22 MPI_Comm comm; 23 24 PetscFunctionBegin; 25 PetscCall(PetscObjectGetComm((PetscObject)pmat, &comm)); 26 PetscCheck(component < blocksize, comm, PETSC_ERR_ARG_INCOMP, "Component %" PetscInt_FMT " should be less than block size %" PetscInt_FMT " ", component, blocksize); 27 PetscCall(MatGetOwnershipRange(pmat, &rstart, &rend)); 28 PetscCheck((rend - rstart) % blocksize == 0, comm, PETSC_ERR_ARG_INCOMP, "Block size %" PetscInt_FMT " is inconsistent for [%" PetscInt_FMT ", %" PetscInt_FMT ") ", blocksize, rstart, rend); 29 PetscCall(ISCreateStride(comm, (rend - rstart) / blocksize, rstart + component, blocksize, &isrow)); 30 PetscCall(MatCreateSubMatrix(pmat, isrow, isrow, reuse, submat)); 31 PetscCall(ISDestroy(&isrow)); 32 PetscFunctionReturn(PETSC_SUCCESS); 33 } 34 35 static PetscErrorCode PCHMGExpandInterpolation_Private(Mat subinterp, Mat *interp, PetscInt blocksize) 36 { 37 PetscInt subrstart, subrend, subrowsize, subcolsize, subcstart, subcend, rowsize, colsize; 38 PetscInt subrow, row, nz, *d_nnz, *o_nnz, i, j, dnz, onz, max_nz, *indices; 39 const PetscInt *idx; 40 const PetscScalar *values; 41 MPI_Comm comm; 42 43 PetscFunctionBegin; 44 PetscCall(PetscObjectGetComm((PetscObject)subinterp, &comm)); 45 PetscCall(MatGetOwnershipRange(subinterp, &subrstart, &subrend)); 46 subrowsize = subrend - subrstart; 47 rowsize = subrowsize * blocksize; 48 PetscCall(PetscCalloc2(rowsize, &d_nnz, rowsize, &o_nnz)); 49 PetscCall(MatGetOwnershipRangeColumn(subinterp, &subcstart, &subcend)); 50 subcolsize = subcend - subcstart; 51 colsize = subcolsize * blocksize; 52 max_nz = 0; 53 for (subrow = subrstart; subrow < subrend; subrow++) { 54 PetscCall(MatGetRow(subinterp, subrow, &nz, &idx, NULL)); 55 if (max_nz < nz) max_nz = nz; 56 dnz = 0; 57 onz = 0; 58 for (i = 0; i < nz; i++) { 59 if (idx[i] >= subcstart && idx[i] < subcend) dnz++; 60 else onz++; 61 } 62 for (i = 0; i < blocksize; i++) { 63 d_nnz[(subrow - subrstart) * blocksize + i] = dnz; 64 o_nnz[(subrow - subrstart) * blocksize + i] = onz; 65 } 66 PetscCall(MatRestoreRow(subinterp, subrow, &nz, &idx, NULL)); 67 } 68 PetscCall(MatCreateAIJ(comm, rowsize, colsize, PETSC_DETERMINE, PETSC_DETERMINE, 0, d_nnz, 0, o_nnz, interp)); 69 PetscCall(MatSetOption(*interp, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE)); 70 PetscCall(MatSetOption(*interp, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 71 PetscCall(MatSetOption(*interp, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 72 PetscCall(MatSetFromOptions(*interp)); 73 74 PetscCall(MatSetUp(*interp)); 75 PetscCall(PetscFree2(d_nnz, o_nnz)); 76 PetscCall(PetscMalloc1(max_nz, &indices)); 77 for (subrow = subrstart; subrow < subrend; subrow++) { 78 PetscCall(MatGetRow(subinterp, subrow, &nz, &idx, &values)); 79 for (i = 0; i < blocksize; i++) { 80 row = subrow * blocksize + i; 81 for (j = 0; j < nz; j++) indices[j] = idx[j] * blocksize + i; 82 PetscCall(MatSetValues(*interp, 1, &row, nz, indices, values, INSERT_VALUES)); 83 } 84 PetscCall(MatRestoreRow(subinterp, subrow, &nz, &idx, &values)); 85 } 86 PetscCall(PetscFree(indices)); 87 PetscCall(MatAssemblyBegin(*interp, MAT_FINAL_ASSEMBLY)); 88 PetscCall(MatAssemblyEnd(*interp, MAT_FINAL_ASSEMBLY)); 89 PetscFunctionReturn(PETSC_SUCCESS); 90 } 91 92 static PetscErrorCode PCSetUp_HMG(PC pc) 93 { 94 Mat PA, submat; 95 PC_MG *mg = (PC_MG *)pc->data; 96 PC_HMG *hmg = (PC_HMG *)mg->innerctx; 97 MPI_Comm comm; 98 PetscInt level; 99 PetscInt num_levels; 100 Mat *operators, *interpolations; 101 PetscInt blocksize; 102 const char *prefix; 103 PCMGGalerkinType galerkin; 104 105 PetscFunctionBegin; 106 PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 107 if (pc->setupcalled) { 108 if (hmg->reuseinterp) { 109 /* If we did not use Galerkin in the last call or we have a different sparsity pattern now, 110 * we have to build from scratch 111 * */ 112 PetscCall(PCMGGetGalerkin(pc, &galerkin)); 113 if (galerkin == PC_MG_GALERKIN_NONE || pc->flag != SAME_NONZERO_PATTERN) pc->setupcalled = PETSC_FALSE; 114 PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_PMAT)); 115 PetscCall(PCSetUp_MG(pc)); 116 PetscFunctionReturn(PETSC_SUCCESS); 117 } else { 118 PetscCall(PCReset_MG(pc)); 119 pc->setupcalled = PETSC_FALSE; 120 } 121 } 122 123 /* Create an inner PC (GAMG or HYPRE) */ 124 if (!hmg->innerpc) { 125 PetscCall(PCCreate(comm, &hmg->innerpc)); 126 /* If users do not set an inner pc type, we need to set a default value */ 127 if (!hmg->innerpctype) { 128 /* If hypre is available, use hypre, otherwise, use gamg */ 129 #if PetscDefined(HAVE_HYPRE) 130 PetscCall(PetscStrallocpy(PCHYPRE, &hmg->innerpctype)); 131 #else 132 PetscCall(PetscStrallocpy(PCGAMG, &hmg->innerpctype)); 133 #endif 134 } 135 PetscCall(PCSetType(hmg->innerpc, hmg->innerpctype)); 136 } 137 PetscCall(PCGetOperators(pc, NULL, &PA)); 138 /* Users need to correctly set a block size of matrix in order to use subspace coarsening */ 139 PetscCall(MatGetBlockSize(PA, &blocksize)); 140 if (blocksize <= 1) hmg->subcoarsening = PETSC_FALSE; 141 /* Extract a submatrix for constructing subinterpolations */ 142 if (hmg->subcoarsening) { 143 PetscCall(PCHMGExtractSubMatrix_Private(PA, &submat, MAT_INITIAL_MATRIX, hmg->component, blocksize)); 144 PA = submat; 145 } 146 PetscCall(PCSetOperators(hmg->innerpc, PA, PA)); 147 if (hmg->subcoarsening) PetscCall(MatDestroy(&PA)); 148 /* Setup inner PC correctly. During this step, matrix will be coarsened */ 149 PetscCall(PCSetUseAmat(hmg->innerpc, PETSC_FALSE)); 150 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 151 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)hmg->innerpc, prefix)); 152 PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)hmg->innerpc, "hmg_inner_")); 153 PetscCall(PCSetFromOptions(hmg->innerpc)); 154 PetscCall(PCSetUp(hmg->innerpc)); 155 156 /* Obtain interpolations IN PLACE. For BoomerAMG, (I,J,data) is reused to avoid memory overhead */ 157 PetscCall(PCGetInterpolations(hmg->innerpc, &num_levels, &interpolations)); 158 /* We can reuse the coarse operators when we do the full space coarsening */ 159 if (!hmg->subcoarsening) PetscCall(PCGetCoarseOperators(hmg->innerpc, &num_levels, &operators)); 160 161 PetscCall(PCDestroy(&hmg->innerpc)); 162 hmg->innerpc = NULL; 163 PetscCall(PCMGSetLevels_MG(pc, num_levels, NULL)); 164 /* Set coarse matrices and interpolations to PCMG */ 165 for (level = num_levels - 1; level > 0; level--) { 166 Mat P = NULL, pmat = NULL; 167 Vec b, x, r; 168 if (hmg->subcoarsening) { 169 if (hmg->usematmaij) { 170 PetscCall(MatCreateMAIJ(interpolations[level - 1], blocksize, &P)); 171 PetscCall(MatDestroy(&interpolations[level - 1])); 172 } else { 173 /* Grow interpolation. In the future, we should use MAIJ */ 174 PetscCall(PCHMGExpandInterpolation_Private(interpolations[level - 1], &P, blocksize)); 175 PetscCall(MatDestroy(&interpolations[level - 1])); 176 } 177 } else { 178 P = interpolations[level - 1]; 179 } 180 PetscCall(MatCreateVecs(P, &b, &r)); 181 PetscCall(PCMGSetInterpolation(pc, level, P)); 182 PetscCall(PCMGSetRestriction(pc, level, P)); 183 PetscCall(MatDestroy(&P)); 184 /* We reuse the matrices when we do not do subspace coarsening */ 185 if ((level - 1) >= 0 && !hmg->subcoarsening) { 186 pmat = operators[level - 1]; 187 PetscCall(PCMGSetOperators(pc, level - 1, pmat, pmat)); 188 PetscCall(MatDestroy(&pmat)); 189 } 190 PetscCall(PCMGSetRhs(pc, level - 1, b)); 191 192 PetscCall(PCMGSetR(pc, level, r)); 193 PetscCall(VecDestroy(&r)); 194 195 PetscCall(VecDuplicate(b, &x)); 196 PetscCall(PCMGSetX(pc, level - 1, x)); 197 PetscCall(VecDestroy(&x)); 198 PetscCall(VecDestroy(&b)); 199 } 200 PetscCall(PetscFree(interpolations)); 201 if (!hmg->subcoarsening) PetscCall(PetscFree(operators)); 202 /* Turn Galerkin off when we already have coarse operators */ 203 PetscCall(PCMGSetGalerkin(pc, hmg->subcoarsening ? PC_MG_GALERKIN_PMAT : PC_MG_GALERKIN_NONE)); 204 PetscCall(PCSetDM(pc, NULL)); 205 PetscCall(PCSetUseAmat(pc, PETSC_FALSE)); 206 PetscObjectOptionsBegin((PetscObject)pc); 207 PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */ 208 PetscOptionsEnd(); 209 PetscCall(PCSetUp_MG(pc)); 210 PetscFunctionReturn(PETSC_SUCCESS); 211 } 212 213 static PetscErrorCode PCDestroy_HMG(PC pc) 214 { 215 PC_MG *mg = (PC_MG *)pc->data; 216 PC_HMG *hmg = (PC_HMG *)mg->innerctx; 217 218 PetscFunctionBegin; 219 PetscCall(PCDestroy(&hmg->innerpc)); 220 PetscCall(PetscFree(hmg->innerpctype)); 221 PetscCall(PetscFree(hmg)); 222 PetscCall(PCDestroy_MG(pc)); 223 224 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetReuseInterpolation_C", NULL)); 225 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetUseSubspaceCoarsening_C", NULL)); 226 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetInnerPCType_C", NULL)); 227 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetCoarseningComponent_C", NULL)); 228 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGUseMatMAIJ_C", NULL)); 229 PetscFunctionReturn(PETSC_SUCCESS); 230 } 231 232 static PetscErrorCode PCView_HMG(PC pc, PetscViewer viewer) 233 { 234 PC_MG *mg = (PC_MG *)pc->data; 235 PC_HMG *hmg = (PC_HMG *)mg->innerctx; 236 PetscBool iascii; 237 238 PetscFunctionBegin; 239 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 240 if (iascii) { 241 PetscCall(PetscViewerASCIIPrintf(viewer, " Reuse interpolation: %s\n", hmg->reuseinterp ? "true" : "false")); 242 PetscCall(PetscViewerASCIIPrintf(viewer, " Use subspace coarsening: %s\n", hmg->subcoarsening ? "true" : "false")); 243 PetscCall(PetscViewerASCIIPrintf(viewer, " Coarsening component: %" PetscInt_FMT " \n", hmg->component)); 244 PetscCall(PetscViewerASCIIPrintf(viewer, " Use MatMAIJ: %s \n", hmg->usematmaij ? "true" : "false")); 245 PetscCall(PetscViewerASCIIPrintf(viewer, " Inner PC type: %s \n", hmg->innerpctype)); 246 } 247 PetscCall(PCView_MG(pc, viewer)); 248 PetscFunctionReturn(PETSC_SUCCESS); 249 } 250 251 static PetscErrorCode PCSetFromOptions_HMG(PC pc, PetscOptionItems PetscOptionsObject) 252 { 253 PC_MG *mg = (PC_MG *)pc->data; 254 PC_HMG *hmg = (PC_HMG *)mg->innerctx; 255 256 PetscFunctionBegin; 257 PetscOptionsHeadBegin(PetscOptionsObject, "HMG"); 258 PetscCall(PetscOptionsBool("-pc_hmg_reuse_interpolation", "Reuse the interpolation operators when possible (cheaper, weaker when matrix entries change a lot)", "PCHMGSetReuseInterpolation", hmg->reuseinterp, &hmg->reuseinterp, NULL)); 259 PetscCall(PetscOptionsBool("-pc_hmg_use_subspace_coarsening", "Use the subspace coarsening to compute the interpolations", "PCHMGSetUseSubspaceCoarsening", hmg->subcoarsening, &hmg->subcoarsening, NULL)); 260 PetscCall(PetscOptionsBool("-pc_hmg_use_matmaij", "Use MatMAIJ store interpolation for saving memory", "PCHMGSetInnerPCType", hmg->usematmaij, &hmg->usematmaij, NULL)); 261 PetscCall(PetscOptionsInt("-pc_hmg_coarsening_component", "Which component is chosen for the subspace-based coarsening algorithm", "PCHMGSetCoarseningComponent", hmg->component, &hmg->component, NULL)); 262 PetscOptionsHeadEnd(); 263 PetscFunctionReturn(PETSC_SUCCESS); 264 } 265 266 static PetscErrorCode PCHMGSetReuseInterpolation_HMG(PC pc, PetscBool reuse) 267 { 268 PC_MG *mg = (PC_MG *)pc->data; 269 PC_HMG *hmg = (PC_HMG *)mg->innerctx; 270 271 PetscFunctionBegin; 272 hmg->reuseinterp = reuse; 273 PetscFunctionReturn(PETSC_SUCCESS); 274 } 275 276 /*@ 277 PCHMGSetReuseInterpolation - Reuse the interpolation matrices in `PCHMG` after changing the matrices numerical values 278 279 Logically Collective 280 281 Input Parameters: 282 + pc - the `PCHMG` context 283 - reuse - `PETSC_TRUE` indicates that `PCHMG` will reuse the interpolations 284 285 Options Database Key: 286 . -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time. 287 288 Level: beginner 289 290 Note: 291 This decreases the set up time of the `PC` significantly but may slow the convergence of the iterative method, `KSP`, that is using the `PCHMG` 292 293 .seealso: [](ch_ksp), `PCHMG`, `PCGAMG`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetCoarseningComponent()`, `PCHMGSetInnerPCType()` 294 @*/ 295 PetscErrorCode PCHMGSetReuseInterpolation(PC pc, PetscBool reuse) 296 { 297 PetscFunctionBegin; 298 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 299 PetscUseMethod(pc, "PCHMGSetReuseInterpolation_C", (PC, PetscBool), (pc, reuse)); 300 PetscFunctionReturn(PETSC_SUCCESS); 301 } 302 303 static PetscErrorCode PCHMGSetUseSubspaceCoarsening_HMG(PC pc, PetscBool subspace) 304 { 305 PC_MG *mg = (PC_MG *)pc->data; 306 PC_HMG *hmg = (PC_HMG *)mg->innerctx; 307 308 PetscFunctionBegin; 309 hmg->subcoarsening = subspace; 310 PetscFunctionReturn(PETSC_SUCCESS); 311 } 312 313 /*@ 314 PCHMGSetUseSubspaceCoarsening - Use subspace coarsening in `PCHMG` 315 316 Logically Collective 317 318 Input Parameters: 319 + pc - the `PCHMG` context 320 - subspace - `PETSC_TRUE` indicates that `PCHMG` will use the subspace coarsening 321 322 Options Database Key: 323 . -pc_hmg_use_subspace_coarsening <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix). 324 325 Level: beginner 326 327 .seealso: [](ch_ksp), `PCHMG`, `PCHMGSetReuseInterpolation()`, `PCHMGSetCoarseningComponent()`, `PCHMGSetInnerPCType()` 328 @*/ 329 PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC pc, PetscBool subspace) 330 { 331 PetscFunctionBegin; 332 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 333 PetscUseMethod(pc, "PCHMGSetUseSubspaceCoarsening_C", (PC, PetscBool), (pc, subspace)); 334 PetscFunctionReturn(PETSC_SUCCESS); 335 } 336 337 static PetscErrorCode PCHMGSetInnerPCType_HMG(PC pc, PCType type) 338 { 339 PC_MG *mg = (PC_MG *)pc->data; 340 PC_HMG *hmg = (PC_HMG *)mg->innerctx; 341 342 PetscFunctionBegin; 343 PetscCall(PetscStrallocpy(type, &hmg->innerpctype)); 344 PetscFunctionReturn(PETSC_SUCCESS); 345 } 346 347 /*@ 348 PCHMGSetInnerPCType - Set an inner `PC` type to be used in the `PCHMG` preconditioner. That is the method used to compute 349 the hierarchy of restriction operators. 350 351 Logically Collective 352 353 Input Parameters: 354 + pc - the `PCHMG` context 355 - type - `PCHYPRE` or `PCGAMG` coarsening algorithm 356 357 Options Database Key: 358 . -hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix 359 360 Level: beginner 361 362 .seealso: [](ch_ksp), `PCHMG`, `PCType`, `PCHMGSetReuseInterpolation()`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetCoarseningComponent()` 363 @*/ 364 PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type) 365 { 366 PetscFunctionBegin; 367 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 368 PetscUseMethod(pc, "PCHMGSetInnerPCType_C", (PC, PCType), (pc, type)); 369 PetscFunctionReturn(PETSC_SUCCESS); 370 } 371 372 static PetscErrorCode PCHMGSetCoarseningComponent_HMG(PC pc, PetscInt component) 373 { 374 PC_MG *mg = (PC_MG *)pc->data; 375 PC_HMG *hmg = (PC_HMG *)mg->innerctx; 376 377 PetscFunctionBegin; 378 hmg->component = component; 379 PetscFunctionReturn(PETSC_SUCCESS); 380 } 381 382 /*@ 383 PCHMGSetCoarseningComponent - Set which component of the PDE is used for the subspace-based coarsening algorithm in the preconditioner `PCHMG` 384 385 Logically Collective 386 387 Input Parameters: 388 + pc - the `PCHMG` context 389 - component - which component `PC` will coarsen 390 391 Options Database Key: 392 . -pc_hmg_coarsening_component <i> - Which component is chosen for the subspace-based coarsening algorithm 393 394 Level: beginner 395 396 Note: 397 By default it uses the first component 398 399 .seealso: [](ch_ksp), `PCHMG`, `PCType`, `PCGAMG`, `PCHMGSetReuseInterpolation()`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetInnerPCType()` 400 @*/ 401 PetscErrorCode PCHMGSetCoarseningComponent(PC pc, PetscInt component) 402 { 403 PetscFunctionBegin; 404 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 405 PetscUseMethod(pc, "PCHMGSetCoarseningComponent_C", (PC, PetscInt), (pc, component)); 406 PetscFunctionReturn(PETSC_SUCCESS); 407 } 408 409 static PetscErrorCode PCHMGUseMatMAIJ_HMG(PC pc, PetscBool usematmaij) 410 { 411 PC_MG *mg = (PC_MG *)pc->data; 412 PC_HMG *hmg = (PC_HMG *)mg->innerctx; 413 414 PetscFunctionBegin; 415 hmg->usematmaij = usematmaij; 416 PetscFunctionReturn(PETSC_SUCCESS); 417 } 418 419 /*@ 420 PCHMGUseMatMAIJ - Set a flag that indicates if or not to use `MATMAIJ` for the interpolation matrices to save memory 421 422 Logically Collective 423 424 Input Parameters: 425 + pc - the `PCHMG` context 426 - usematmaij - `PETSC_TRUE` (default) to use `MATMAIJ` for interpolations. 427 428 Options Database Key: 429 . -pc_hmg_use_matmaij - <true | false > 430 431 Level: beginner 432 433 .seealso: [](ch_ksp), `PCHMG`, `PCType`, `PCGAMG` 434 @*/ 435 PetscErrorCode PCHMGUseMatMAIJ(PC pc, PetscBool usematmaij) 436 { 437 PetscFunctionBegin; 438 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 439 PetscUseMethod(pc, "PCHMGUseMatMAIJ_C", (PC, PetscBool), (pc, usematmaij)); 440 PetscFunctionReturn(PETSC_SUCCESS); 441 } 442 443 /*MC 444 PCHMG - Preconditioner for multiple component PDE problems that constructs a hierarchy of restriction operators to coarse grid problems using the submatrix of 445 a single component with either `PCHYPRE` or `PCGAMG`. The same restriction operators are then used for each of the components of the PDE within the `PCMG` 446 multigrid preconditioner. This results in a much more efficient to build and apply preconditioner than using `PCGAMG` on the entire system {cite}`kong2020highly`. 447 448 Options Database Keys: 449 + -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations for new matrix values or rebuild the interpolation. This can save compute time. 450 . -pc_hmg_use_subspace_coarsening <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix, or coarsen on the full matrix). 451 . -hmg_inner_pc_type <hypre, gamg, ...> - What method to use to generate the hierarchy of restriction operators 452 - -pc_hmg_use_matmaij <true | false> - Whether or not to use `MATMAIJ` for multicomponent problems for saving memory 453 454 Level: intermediate 455 456 Note: 457 `MatSetBlockSize()` must be called on the linear system matrix to set the number of components of the PDE. 458 459 .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMG`, `PCHYPRE`, `PCHMG`, `PCGetCoarseOperators()`, `PCGetInterpolations()`, 460 `PCHMGSetReuseInterpolation()`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetInnerPCType()`, `PCGAMG` 461 M*/ 462 PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc) 463 { 464 PC_HMG *hmg; 465 PC_MG *mg; 466 467 PetscFunctionBegin; 468 /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */ 469 PetscTryTypeMethod(pc, destroy); 470 pc->data = NULL; 471 PetscCall(PetscFree(((PetscObject)pc)->type_name)); 472 473 PetscCall(PCSetType(pc, PCMG)); 474 PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCHMG)); 475 PetscCall(PetscNew(&hmg)); 476 477 mg = (PC_MG *)pc->data; 478 mg->innerctx = hmg; 479 hmg->reuseinterp = PETSC_FALSE; 480 hmg->subcoarsening = PETSC_FALSE; 481 hmg->usematmaij = PETSC_TRUE; 482 hmg->component = 0; 483 hmg->innerpc = NULL; 484 485 pc->ops->setfromoptions = PCSetFromOptions_HMG; 486 pc->ops->view = PCView_HMG; 487 pc->ops->destroy = PCDestroy_HMG; 488 pc->ops->setup = PCSetUp_HMG; 489 490 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetReuseInterpolation_C", PCHMGSetReuseInterpolation_HMG)); 491 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetUseSubspaceCoarsening_C", PCHMGSetUseSubspaceCoarsening_HMG)); 492 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetInnerPCType_C", PCHMGSetInnerPCType_HMG)); 493 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetCoarseningComponent_C", PCHMGSetCoarseningComponent_HMG)); 494 PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGUseMatMAIJ_C", PCHMGUseMatMAIJ_HMG)); 495 PetscFunctionReturn(PETSC_SUCCESS); 496 } 497