1360ee056SFande Kong #include <petscdm.h> 2eec179cfSJacob Faibussowitsch #include <petsc/private/hashmapi.h> 3360ee056SFande Kong #include <petsc/private/matimpl.h> 4360ee056SFande Kong #include <petsc/private/pcmgimpl.h> 5360ee056SFande Kong #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 6360ee056SFande Kong 7360ee056SFande Kong typedef struct { 88a2c336bSFande Kong PC innerpc; /* A MG inner PC (Hypre or PCGAMG) to setup interpolations and coarse operators */ 98a2c336bSFande Kong char *innerpctype; /* PCGAMG or PCHYPRE */ 108a2c336bSFande Kong PetscBool reuseinterp; /* A flag indicates if or not to reuse the interpolations */ 1149c604d5SFande Kong PetscBool subcoarsening; /* If or not to use a subspace-based coarsening algorithm */ 1249c604d5SFande Kong PetscBool usematmaij; /* If or not to use MatMAIJ for saving memory */ 1349c604d5SFande Kong PetscInt component; /* Which subspace is used for the subspace-based coarsening algorithm? */ 14360ee056SFande Kong } PC_HMG; 15360ee056SFande Kong 1666976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_HMG(PC, PetscOptionItems *); 17360ee056SFande Kong PetscErrorCode PCReset_MG(PC); 188a2c336bSFande Kong 19d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGExtractSubMatrix_Private(Mat pmat, Mat *submat, MatReuse reuse, PetscInt component, PetscInt blocksize) 20d71ae5a4SJacob Faibussowitsch { 21360ee056SFande Kong IS isrow; 228a2c336bSFande Kong PetscInt rstart, rend; 23360ee056SFande Kong MPI_Comm comm; 24360ee056SFande Kong 25360ee056SFande Kong PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pmat, &comm)); 2763a3b9bcSJacob Faibussowitsch PetscCheck(component < blocksize, comm, PETSC_ERR_ARG_INCOMP, "Component %" PetscInt_FMT " should be less than block size %" PetscInt_FMT " ", component, blocksize); 289566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pmat, &rstart, &rend)); 2963a3b9bcSJacob Faibussowitsch PetscCheck((rend - rstart) % blocksize == 0, comm, PETSC_ERR_ARG_INCOMP, "Block size %" PetscInt_FMT " is inconsistent for [%" PetscInt_FMT ", %" PetscInt_FMT ") ", blocksize, rstart, rend); 309566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, (rend - rstart) / blocksize, rstart + component, blocksize, &isrow)); 319566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pmat, isrow, isrow, reuse, submat)); 329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34360ee056SFande Kong } 35360ee056SFande Kong 36d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGExpandInterpolation_Private(Mat subinterp, Mat *interp, PetscInt blocksize) 37d71ae5a4SJacob Faibussowitsch { 38360ee056SFande Kong PetscInt subrstart, subrend, subrowsize, subcolsize, subcstart, subcend, rowsize, colsize; 398a2c336bSFande Kong PetscInt subrow, row, nz, *d_nnz, *o_nnz, i, j, dnz, onz, max_nz, *indices; 408a2c336bSFande Kong const PetscInt *idx; 418a2c336bSFande Kong const PetscScalar *values; 428a2c336bSFande Kong MPI_Comm comm; 43360ee056SFande Kong 44360ee056SFande Kong PetscFunctionBegin; 459566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)subinterp, &comm)); 469566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(subinterp, &subrstart, &subrend)); 478a2c336bSFande Kong subrowsize = subrend - subrstart; 48360ee056SFande Kong rowsize = subrowsize * blocksize; 499566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(rowsize, &d_nnz, rowsize, &o_nnz)); 509566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(subinterp, &subcstart, &subcend)); 518a2c336bSFande Kong subcolsize = subcend - subcstart; 528a2c336bSFande Kong colsize = subcolsize * blocksize; 53360ee056SFande Kong max_nz = 0; 548a2c336bSFande Kong for (subrow = subrstart; subrow < subrend; subrow++) { 559566063dSJacob Faibussowitsch PetscCall(MatGetRow(subinterp, subrow, &nz, &idx, NULL)); 56360ee056SFande Kong if (max_nz < nz) max_nz = nz; 579371c9d4SSatish Balay dnz = 0; 589371c9d4SSatish Balay onz = 0; 59360ee056SFande Kong for (i = 0; i < nz; i++) { 608a2c336bSFande Kong if (idx[i] >= subcstart && idx[i] < subcend) dnz++; 618a2c336bSFande Kong else onz++; 62360ee056SFande Kong } 63360ee056SFande Kong for (i = 0; i < blocksize; i++) { 64360ee056SFande Kong d_nnz[(subrow - subrstart) * blocksize + i] = dnz; 65360ee056SFande Kong o_nnz[(subrow - subrstart) * blocksize + i] = onz; 66360ee056SFande Kong } 679566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(subinterp, subrow, &nz, &idx, NULL)); 68360ee056SFande Kong } 699566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(comm, rowsize, colsize, PETSC_DETERMINE, PETSC_DETERMINE, 0, d_nnz, 0, o_nnz, interp)); 709566063dSJacob Faibussowitsch PetscCall(MatSetOption(*interp, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE)); 719566063dSJacob Faibussowitsch PetscCall(MatSetOption(*interp, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE)); 729566063dSJacob Faibussowitsch PetscCall(MatSetOption(*interp, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 739566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(*interp)); 74360ee056SFande Kong 759566063dSJacob Faibussowitsch PetscCall(MatSetUp(*interp)); 769566063dSJacob Faibussowitsch PetscCall(PetscFree2(d_nnz, o_nnz)); 779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(max_nz, &indices)); 788a2c336bSFande Kong for (subrow = subrstart; subrow < subrend; subrow++) { 799566063dSJacob Faibussowitsch PetscCall(MatGetRow(subinterp, subrow, &nz, &idx, &values)); 80360ee056SFande Kong for (i = 0; i < blocksize; i++) { 81360ee056SFande Kong row = subrow * blocksize + i; 82ad540459SPierre Jolivet for (j = 0; j < nz; j++) indices[j] = idx[j] * blocksize + i; 839566063dSJacob Faibussowitsch PetscCall(MatSetValues(*interp, 1, &row, nz, indices, values, INSERT_VALUES)); 84360ee056SFande Kong } 859566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(subinterp, subrow, &nz, &idx, &values)); 86360ee056SFande Kong } 879566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 889566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*interp, MAT_FINAL_ASSEMBLY)); 899566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*interp, MAT_FINAL_ASSEMBLY)); 903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91360ee056SFande Kong } 92360ee056SFande Kong 9366976f2fSJacob Faibussowitsch static PetscErrorCode PCSetUp_HMG(PC pc) 94d71ae5a4SJacob Faibussowitsch { 95360ee056SFande Kong Mat PA, submat; 96360ee056SFande Kong PC_MG *mg = (PC_MG *)pc->data; 97360ee056SFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx; 98360ee056SFande Kong MPI_Comm comm; 99360ee056SFande Kong PetscInt level; 100360ee056SFande Kong PetscInt num_levels; 1018a2c336bSFande Kong Mat *operators, *interpolations; 1028a2c336bSFande Kong PetscInt blocksize; 103fd2dd295SFande Kong const char *prefix; 10407a4832bSFande Kong PCMGGalerkinType galerkin; 105360ee056SFande Kong 106360ee056SFande Kong PetscFunctionBegin; 1079566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 108360ee056SFande Kong if (pc->setupcalled) { 1095f36d8f2SFande Kong if (hmg->reuseinterp) { 1105f36d8f2SFande Kong /* If we did not use Galerkin in the last call or we have a different sparsity pattern now, 1115f36d8f2SFande Kong * we have to build from scratch 11207a4832bSFande Kong * */ 1139566063dSJacob Faibussowitsch PetscCall(PCMGGetGalerkin(pc, &galerkin)); 1145f36d8f2SFande Kong if (galerkin == PC_MG_GALERKIN_NONE || pc->flag != SAME_NONZERO_PATTERN) pc->setupcalled = PETSC_FALSE; 1159566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_PMAT)); 1169566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 118360ee056SFande Kong } else { 1199566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc)); 120360ee056SFande Kong pc->setupcalled = PETSC_FALSE; 121360ee056SFande Kong } 122360ee056SFande Kong } 123360ee056SFande Kong 1248a2c336bSFande Kong /* Create an inner PC (GAMG or HYPRE) */ 1258a2c336bSFande Kong if (!hmg->innerpc) { 1269566063dSJacob Faibussowitsch PetscCall(PCCreate(comm, &hmg->innerpc)); 127fd2dd295SFande Kong /* If users do not set an inner pc type, we need to set a default value */ 128fd2dd295SFande Kong if (!hmg->innerpctype) { 129fd2dd295SFande Kong /* If hypre is available, use hypre, otherwise, use gamg */ 130a51cd166SPierre Jolivet #if PetscDefined(HAVE_HYPRE) 131*f4f49eeaSPierre Jolivet PetscCall(PetscStrallocpy(PCHYPRE, &hmg->innerpctype)); 132fd2dd295SFande Kong #else 133*f4f49eeaSPierre Jolivet PetscCall(PetscStrallocpy(PCGAMG, &hmg->innerpctype)); 134fd2dd295SFande Kong #endif 1358a2c336bSFande Kong } 1369566063dSJacob Faibussowitsch PetscCall(PCSetType(hmg->innerpc, hmg->innerpctype)); 137360ee056SFande Kong } 1389566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &PA)); 1398a2c336bSFande Kong /* Users need to correctly set a block size of matrix in order to use subspace coarsening */ 1409566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(PA, &blocksize)); 1418a2c336bSFande Kong if (blocksize <= 1) hmg->subcoarsening = PETSC_FALSE; 1428a2c336bSFande Kong /* Extract a submatrix for constructing subinterpolations */ 1438a2c336bSFande Kong if (hmg->subcoarsening) { 1449566063dSJacob Faibussowitsch PetscCall(PCHMGExtractSubMatrix_Private(PA, &submat, MAT_INITIAL_MATRIX, hmg->component, blocksize)); 145360ee056SFande Kong PA = submat; 146360ee056SFande Kong } 1479566063dSJacob Faibussowitsch PetscCall(PCSetOperators(hmg->innerpc, PA, PA)); 14848a46eb9SPierre Jolivet if (hmg->subcoarsening) PetscCall(MatDestroy(&PA)); 1498a2c336bSFande Kong /* Setup inner PC correctly. During this step, matrix will be coarsened */ 1509566063dSJacob Faibussowitsch PetscCall(PCSetUseAmat(hmg->innerpc, PETSC_FALSE)); 1519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix)); 1529566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)hmg->innerpc, prefix)); 1539566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)hmg->innerpc, "hmg_inner_")); 1549566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(hmg->innerpc)); 1559566063dSJacob Faibussowitsch PetscCall(PCSetUp(hmg->innerpc)); 156360ee056SFande Kong 1578a2c336bSFande Kong /* Obtain interpolations IN PLACE. For BoomerAMG, (I,J,data) is reused to avoid memory overhead */ 1589566063dSJacob Faibussowitsch PetscCall(PCGetInterpolations(hmg->innerpc, &num_levels, &interpolations)); 1598a2c336bSFande Kong /* We can reuse the coarse operators when we do the full space coarsening */ 16048a46eb9SPierre Jolivet if (!hmg->subcoarsening) PetscCall(PCGetCoarseOperators(hmg->innerpc, &num_levels, &operators)); 161360ee056SFande Kong 1629566063dSJacob Faibussowitsch PetscCall(PCDestroy(&hmg->innerpc)); 1630a545947SLisandro Dalcin hmg->innerpc = NULL; 1649566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels_MG(pc, num_levels, NULL)); 1658a2c336bSFande Kong /* Set coarse matrices and interpolations to PCMG */ 166360ee056SFande Kong for (level = num_levels - 1; level > 0; level--) { 1670a545947SLisandro Dalcin Mat P = NULL, pmat = NULL; 168360ee056SFande Kong Vec b, x, r; 1698a2c336bSFande Kong if (hmg->subcoarsening) { 1704bb91820SFande Kong if (hmg->usematmaij) { 1719566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(interpolations[level - 1], blocksize, &P)); 1729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&interpolations[level - 1])); 1734bb91820SFande Kong } else { 1748a2c336bSFande Kong /* Grow interpolation. In the future, we should use MAIJ */ 1759566063dSJacob Faibussowitsch PetscCall(PCHMGExpandInterpolation_Private(interpolations[level - 1], &P, blocksize)); 1769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&interpolations[level - 1])); 1774bb91820SFande Kong } 178360ee056SFande Kong } else { 1798a2c336bSFande Kong P = interpolations[level - 1]; 180360ee056SFande Kong } 1819566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(P, &b, &r)); 1829566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, level, P)); 1839566063dSJacob Faibussowitsch PetscCall(PCMGSetRestriction(pc, level, P)); 1849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P)); 1858a2c336bSFande Kong /* We reuse the matrices when we do not do subspace coarsening */ 1868a2c336bSFande Kong if ((level - 1) >= 0 && !hmg->subcoarsening) { 1878a2c336bSFande Kong pmat = operators[level - 1]; 1889566063dSJacob Faibussowitsch PetscCall(PCMGSetOperators(pc, level - 1, pmat, pmat)); 1899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pmat)); 190360ee056SFande Kong } 1919566063dSJacob Faibussowitsch PetscCall(PCMGSetRhs(pc, level - 1, b)); 192360ee056SFande Kong 1939566063dSJacob Faibussowitsch PetscCall(PCMGSetR(pc, level, r)); 1949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r)); 195360ee056SFande Kong 1969566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b, &x)); 1979566063dSJacob Faibussowitsch PetscCall(PCMGSetX(pc, level - 1, x)); 1989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 200360ee056SFande Kong } 2019566063dSJacob Faibussowitsch PetscCall(PetscFree(interpolations)); 20248a46eb9SPierre Jolivet if (!hmg->subcoarsening) PetscCall(PetscFree(operators)); 2038a2c336bSFande Kong /* Turn Galerkin off when we already have coarse operators */ 2049566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, hmg->subcoarsening ? PC_MG_GALERKIN_PMAT : PC_MG_GALERKIN_NONE)); 2059566063dSJacob Faibussowitsch PetscCall(PCSetDM(pc, NULL)); 2069566063dSJacob Faibussowitsch PetscCall(PCSetUseAmat(pc, PETSC_FALSE)); 207d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)pc); 208dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */ 209d0609cedSBarry Smith PetscOptionsEnd(); 2109566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc)); 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 212360ee056SFande Kong } 213360ee056SFande Kong 21466976f2fSJacob Faibussowitsch static PetscErrorCode PCDestroy_HMG(PC pc) 215d71ae5a4SJacob Faibussowitsch { 216360ee056SFande Kong PC_MG *mg = (PC_MG *)pc->data; 2178a2c336bSFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx; 218360ee056SFande Kong 219360ee056SFande Kong PetscFunctionBegin; 2209566063dSJacob Faibussowitsch PetscCall(PCDestroy(&hmg->innerpc)); 2219566063dSJacob Faibussowitsch PetscCall(PetscFree(hmg->innerpctype)); 2229566063dSJacob Faibussowitsch PetscCall(PetscFree(hmg)); 2239566063dSJacob Faibussowitsch PetscCall(PCDestroy_MG(pc)); 224fd2dd295SFande Kong 2259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetReuseInterpolation_C", NULL)); 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetUseSubspaceCoarsening_C", NULL)); 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetInnerPCType_C", NULL)); 2289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetCoarseningComponent_C", NULL)); 2292e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGUseMatMAIJ_C", NULL)); 2303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 231360ee056SFande Kong } 232360ee056SFande Kong 23366976f2fSJacob Faibussowitsch static PetscErrorCode PCView_HMG(PC pc, PetscViewer viewer) 234d71ae5a4SJacob Faibussowitsch { 235360ee056SFande Kong PC_MG *mg = (PC_MG *)pc->data; 2368a2c336bSFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx; 237360ee056SFande Kong PetscBool iascii; 238360ee056SFande Kong 239360ee056SFande Kong PetscFunctionBegin; 2409566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 241360ee056SFande Kong if (iascii) { 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Reuse interpolation: %s\n", hmg->reuseinterp ? "true" : "false")); 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Use subspace coarsening: %s\n", hmg->subcoarsening ? "true" : "false")); 24463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Coarsening component: %" PetscInt_FMT " \n", hmg->component)); 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Use MatMAIJ: %s \n", hmg->usematmaij ? "true" : "false")); 2469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Inner PC type: %s \n", hmg->innerpctype)); 247360ee056SFande Kong } 2489566063dSJacob Faibussowitsch PetscCall(PCView_MG(pc, viewer)); 2493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 250360ee056SFande Kong } 251360ee056SFande Kong 25266976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_HMG(PC pc, PetscOptionItems *PetscOptionsObject) 253d71ae5a4SJacob Faibussowitsch { 254360ee056SFande Kong PC_MG *mg = (PC_MG *)pc->data; 2558a2c336bSFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx; 256360ee056SFande Kong 257360ee056SFande Kong PetscFunctionBegin; 258d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "HMG"); 2599566063dSJacob Faibussowitsch 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)); 2609566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_hmg_use_subspace_coarsening", "Use the subspace coarsening to compute the interpolations", "PCHMGSetUseSubspaceCoarsening", hmg->subcoarsening, &hmg->subcoarsening, NULL)); 2619566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_hmg_use_matmaij", "Use MatMAIJ store interpolation for saving memory", "PCHMGSetInnerPCType", hmg->usematmaij, &hmg->usematmaij, NULL)); 2629566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_hmg_coarsening_component", "Which component is chosen for the subspace-based coarsening algorithm", "PCHMGSetCoarseningComponent", hmg->component, &hmg->component, NULL)); 263d0609cedSBarry Smith PetscOptionsHeadEnd(); 2643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 265360ee056SFande Kong } 266360ee056SFande Kong 267d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGSetReuseInterpolation_HMG(PC pc, PetscBool reuse) 268d71ae5a4SJacob Faibussowitsch { 26907a4832bSFande Kong PC_MG *mg = (PC_MG *)pc->data; 27007a4832bSFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx; 271fd2dd295SFande Kong 272fd2dd295SFande Kong PetscFunctionBegin; 273fd2dd295SFande Kong hmg->reuseinterp = reuse; 2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 275fd2dd295SFande Kong } 276fd2dd295SFande Kong 277b155ba7fSFande Kong /*@ 278f1580f4eSBarry Smith PCHMGSetReuseInterpolation - Reuse the interpolation matrices in `PCHMG` after changing the matrices numerical values 2798a2c336bSFande Kong 280c3339decSBarry Smith Logically Collective 2818a2c336bSFande Kong 2828a2c336bSFande Kong Input Parameters: 283f1580f4eSBarry Smith + pc - the `PCHMG` context 284f1580f4eSBarry Smith - reuse - `PETSC_TRUE` indicates that `PCHMG` will reuse the interpolations 2858a2c336bSFande Kong 286f1580f4eSBarry Smith Options Database Key: 28749c604d5SFande Kong . -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time. 2888a2c336bSFande Kong 2898a2c336bSFande Kong Level: beginner 2908a2c336bSFande Kong 291562efe2eSBarry Smith .seealso: [](ch_ksp), `PCHMG`, `PCGAMG`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetCoarseningComponent()`, `PCHMGSetInnerPCType()` 292b155ba7fSFande Kong @*/ 293d71ae5a4SJacob Faibussowitsch PetscErrorCode PCHMGSetReuseInterpolation(PC pc, PetscBool reuse) 294d71ae5a4SJacob Faibussowitsch { 295fd2dd295SFande Kong PetscFunctionBegin; 296fd2dd295SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 297cac4c232SBarry Smith PetscUseMethod(pc, "PCHMGSetReuseInterpolation_C", (PC, PetscBool), (pc, reuse)); 2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 299fd2dd295SFande Kong } 300fd2dd295SFande Kong 301d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGSetUseSubspaceCoarsening_HMG(PC pc, PetscBool subspace) 302d71ae5a4SJacob Faibussowitsch { 30307a4832bSFande Kong PC_MG *mg = (PC_MG *)pc->data; 30407a4832bSFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx; 3058a2c336bSFande Kong 3068a2c336bSFande Kong PetscFunctionBegin; 307fd2dd295SFande Kong hmg->subcoarsening = subspace; 3083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3098a2c336bSFande Kong } 3108a2c336bSFande Kong 311b155ba7fSFande Kong /*@ 312f1580f4eSBarry Smith PCHMGSetUseSubspaceCoarsening - Use subspace coarsening in `PCHMG` 3138a2c336bSFande Kong 314c3339decSBarry Smith Logically Collective 3158a2c336bSFande Kong 3168a2c336bSFande Kong Input Parameters: 317f1580f4eSBarry Smith + pc - the `PCHMG` context 318feefa0e1SJacob Faibussowitsch - subspace - `PETSC_TRUE` indicates that `PCHMG` will use the subspace coarsening 3198a2c336bSFande Kong 320f1580f4eSBarry Smith Options Database Key: 32149c604d5SFande Kong . -pc_hmg_use_subspace_coarsening <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix). 3228a2c336bSFande Kong 3238a2c336bSFande Kong Level: beginner 3248a2c336bSFande Kong 325562efe2eSBarry Smith .seealso: [](ch_ksp), `PCHMG`, `PCHMGSetReuseInterpolation()`, `PCHMGSetCoarseningComponent()`, `PCHMGSetInnerPCType()` 326b155ba7fSFande Kong @*/ 327d71ae5a4SJacob Faibussowitsch PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC pc, PetscBool subspace) 328d71ae5a4SJacob Faibussowitsch { 3298a2c336bSFande Kong PetscFunctionBegin; 3308a2c336bSFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 331cac4c232SBarry Smith PetscUseMethod(pc, "PCHMGSetUseSubspaceCoarsening_C", (PC, PetscBool), (pc, subspace)); 3323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 333fd2dd295SFande Kong } 334fd2dd295SFande Kong 335d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGSetInnerPCType_HMG(PC pc, PCType type) 336d71ae5a4SJacob Faibussowitsch { 33707a4832bSFande Kong PC_MG *mg = (PC_MG *)pc->data; 33807a4832bSFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx; 339fd2dd295SFande Kong 340fd2dd295SFande Kong PetscFunctionBegin; 341*f4f49eeaSPierre Jolivet PetscCall(PetscStrallocpy(type, &hmg->innerpctype)); 3423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3438a2c336bSFande Kong } 3448a2c336bSFande Kong 345b155ba7fSFande Kong /*@C 346f1580f4eSBarry Smith PCHMGSetInnerPCType - Set an inner `PC` type 3478a2c336bSFande Kong 348c3339decSBarry Smith Logically Collective 3498a2c336bSFande Kong 3508a2c336bSFande Kong Input Parameters: 351f1580f4eSBarry Smith + pc - the `PCHMG` context 352f1580f4eSBarry Smith - type - `PCHYPRE` or `PCGAMG` coarsening algorithm 3538a2c336bSFande Kong 354f1580f4eSBarry Smith Options Database Key: 35549c604d5SFande Kong . -hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix 3568a2c336bSFande Kong 3578a2c336bSFande Kong Level: beginner 3588a2c336bSFande Kong 359562efe2eSBarry Smith .seealso: [](ch_ksp), `PCHMG`, `PCType`, `PCHMGSetReuseInterpolation()`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetCoarseningComponent()` 360b155ba7fSFande Kong @*/ 361d71ae5a4SJacob Faibussowitsch PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type) 362d71ae5a4SJacob Faibussowitsch { 3638a2c336bSFande Kong PetscFunctionBegin; 3648a2c336bSFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 365cac4c232SBarry Smith PetscUseMethod(pc, "PCHMGSetInnerPCType_C", (PC, PCType), (pc, type)); 3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3678a2c336bSFande Kong } 3688a2c336bSFande Kong 369d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGSetCoarseningComponent_HMG(PC pc, PetscInt component) 370d71ae5a4SJacob Faibussowitsch { 37149c604d5SFande Kong PC_MG *mg = (PC_MG *)pc->data; 37249c604d5SFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx; 37349c604d5SFande Kong 37449c604d5SFande Kong PetscFunctionBegin; 37549c604d5SFande Kong hmg->component = component; 3763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 37749c604d5SFande Kong } 37849c604d5SFande Kong 379b155ba7fSFande Kong /*@ 380f1580f4eSBarry Smith PCHMGSetCoarseningComponent - Set which component of the PDE is used for the subspace-based coarsening algorithm 38149c604d5SFande Kong 382c3339decSBarry Smith Logically Collective 38349c604d5SFande Kong 38449c604d5SFande Kong Input Parameters: 385f1580f4eSBarry Smith + pc - the `PCHMG` context 386f1580f4eSBarry Smith - component - which component `PC` will coarsen 38749c604d5SFande Kong 388f1580f4eSBarry Smith Options Database Key: 389f1580f4eSBarry Smith . -pc_hmg_coarsening_component <i> - Which component is chosen for the subspace-based coarsening algorithm 39049c604d5SFande Kong 39149c604d5SFande Kong Level: beginner 39249c604d5SFande Kong 393562efe2eSBarry Smith .seealso: [](ch_ksp), `PCHMG`, `PCType`, `PCGAMG`, `PCHMGSetReuseInterpolation()`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetInnerPCType()` 394b155ba7fSFande Kong @*/ 395d71ae5a4SJacob Faibussowitsch PetscErrorCode PCHMGSetCoarseningComponent(PC pc, PetscInt component) 396d71ae5a4SJacob Faibussowitsch { 39749c604d5SFande Kong PetscFunctionBegin; 39849c604d5SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 399cac4c232SBarry Smith PetscUseMethod(pc, "PCHMGSetCoarseningComponent_C", (PC, PetscInt), (pc, component)); 4003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40149c604d5SFande Kong } 40249c604d5SFande Kong 403d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGUseMatMAIJ_HMG(PC pc, PetscBool usematmaij) 404d71ae5a4SJacob Faibussowitsch { 40549c604d5SFande Kong PC_MG *mg = (PC_MG *)pc->data; 40649c604d5SFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx; 40749c604d5SFande Kong 40849c604d5SFande Kong PetscFunctionBegin; 40949c604d5SFande Kong hmg->usematmaij = usematmaij; 4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41149c604d5SFande Kong } 41249c604d5SFande Kong 413b155ba7fSFande Kong /*@ 414f1580f4eSBarry Smith PCHMGUseMatMAIJ - Set a flag that indicates if or not to use `MATMAIJ` for the interpolation matrices for saving memory 41549c604d5SFande Kong 416c3339decSBarry Smith Logically Collective 41749c604d5SFande Kong 41849c604d5SFande Kong Input Parameters: 419f1580f4eSBarry Smith + pc - the `PCHMG` context 420f1580f4eSBarry Smith - usematmaij - `PETSC_TRUE` (default) to use `MATMAIJ` for interpolations. 42149c604d5SFande Kong 422f1580f4eSBarry Smith Options Database Key: 42349c604d5SFande Kong . -pc_hmg_use_matmaij - <true | false > 42449c604d5SFande Kong 42549c604d5SFande Kong Level: beginner 42649c604d5SFande Kong 427562efe2eSBarry Smith .seealso: [](ch_ksp), `PCHMG`, `PCType`, `PCGAMG` 428b155ba7fSFande Kong @*/ 429d71ae5a4SJacob Faibussowitsch PetscErrorCode PCHMGUseMatMAIJ(PC pc, PetscBool usematmaij) 430d71ae5a4SJacob Faibussowitsch { 43149c604d5SFande Kong PetscFunctionBegin; 43249c604d5SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 433cac4c232SBarry Smith PetscUseMethod(pc, "PCHMGUseMatMAIJ_C", (PC, PetscBool), (pc, usematmaij)); 4343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43549c604d5SFande Kong } 43649c604d5SFande Kong 4378a2c336bSFande Kong /*MC 438f1580f4eSBarry Smith PCHMG - For multiple component PDE problems constructs a hierarchy of restriction operators to coarse grid problems using the submatrix of 439f1580f4eSBarry Smith a single component with either `PCHYPRE` or `PCGAMG`. The same restriction operators are used for each of the components of the PDE with `PCMG` 4401d27aa22SBarry Smith resulting in a much more efficient to build and apply preconditioner than using `PCGAMG` on the entire system {cite}`kong2020highly`. 441360ee056SFande Kong 442360ee056SFande Kong Options Database Keys: 443f1580f4eSBarry Smith + -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations for new matrix values. It can potentially save compute time. 4448a2c336bSFande Kong . -pc_hmg_use_subspace_coarsening <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix). 4451d27aa22SBarry Smith . -hmg_inner_pc_type <hypre, gamg, ...> - What method is used to solve the coarsen matrix 446f1580f4eSBarry Smith - -pc_hmg_use_matmaij <true | false> - Whether or not to use `MATMAIJ` for multicomponent problems for saving memory 447360ee056SFande Kong 448f1580f4eSBarry Smith Level: intermediate 449360ee056SFande Kong 450f1580f4eSBarry Smith Note: 451f1580f4eSBarry Smith `MatSetBlockSize()` must be called on the linear system matrix to set the number of components of the PDE. 452360ee056SFande Kong 453562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMG`, `PCHYPRE`, `PCHMG`, `PCGetCoarseOperators()`, `PCGetInterpolations()`, 454f1580f4eSBarry Smith `PCHMGSetReuseInterpolation()`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetInnerPCType()` 455360ee056SFande Kong M*/ 456d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc) 457d71ae5a4SJacob Faibussowitsch { 458360ee056SFande Kong PC_HMG *hmg; 459360ee056SFande Kong PC_MG *mg; 460360ee056SFande Kong 461360ee056SFande Kong PetscFunctionBegin; 462360ee056SFande Kong /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */ 463dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, destroy); 4640a545947SLisandro Dalcin pc->data = NULL; 4659566063dSJacob Faibussowitsch PetscCall(PetscFree(((PetscObject)pc)->type_name)); 466360ee056SFande Kong 4679566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCMG)); 4689566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCHMG)); 4699566063dSJacob Faibussowitsch PetscCall(PetscNew(&hmg)); 470360ee056SFande Kong 471360ee056SFande Kong mg = (PC_MG *)pc->data; 472360ee056SFande Kong mg->innerctx = hmg; 473360ee056SFande Kong hmg->reuseinterp = PETSC_FALSE; 4748a2c336bSFande Kong hmg->subcoarsening = PETSC_FALSE; 4754bb91820SFande Kong hmg->usematmaij = PETSC_TRUE; 47649c604d5SFande Kong hmg->component = 0; 4778a2c336bSFande Kong hmg->innerpc = NULL; 478360ee056SFande Kong 479360ee056SFande Kong pc->ops->setfromoptions = PCSetFromOptions_HMG; 480360ee056SFande Kong pc->ops->view = PCView_HMG; 481360ee056SFande Kong pc->ops->destroy = PCDestroy_HMG; 482360ee056SFande Kong pc->ops->setup = PCSetUp_HMG; 483fd2dd295SFande Kong 4849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetReuseInterpolation_C", PCHMGSetReuseInterpolation_HMG)); 4859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetUseSubspaceCoarsening_C", PCHMGSetUseSubspaceCoarsening_HMG)); 4869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetInnerPCType_C", PCHMGSetInnerPCType_HMG)); 4879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetCoarseningComponent_C", PCHMGSetCoarseningComponent_HMG)); 4889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGUseMatMAIJ_C", PCHMGUseMatMAIJ_HMG)); 4893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 490360ee056SFande Kong } 491