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
PCHMGExtractSubMatrix_Private(Mat pmat,Mat * submat,MatReuse reuse,PetscInt component,PetscInt blocksize)16d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGExtractSubMatrix_Private(Mat pmat, Mat *submat, MatReuse reuse, PetscInt component, PetscInt blocksize)
17d71ae5a4SJacob Faibussowitsch {
18360ee056SFande Kong IS isrow;
198a2c336bSFande Kong PetscInt rstart, rend;
20360ee056SFande Kong MPI_Comm comm;
21360ee056SFande Kong
22360ee056SFande Kong PetscFunctionBegin;
239566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pmat, &comm));
2463a3b9bcSJacob Faibussowitsch PetscCheck(component < blocksize, comm, PETSC_ERR_ARG_INCOMP, "Component %" PetscInt_FMT " should be less than block size %" PetscInt_FMT " ", component, blocksize);
259566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pmat, &rstart, &rend));
2663a3b9bcSJacob 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);
279566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, (rend - rstart) / blocksize, rstart + component, blocksize, &isrow));
289566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pmat, isrow, isrow, reuse, submat));
299566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow));
303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
31360ee056SFande Kong }
32360ee056SFande Kong
PCHMGExpandInterpolation_Private(Mat subinterp,Mat * interp,PetscInt blocksize)33d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGExpandInterpolation_Private(Mat subinterp, Mat *interp, PetscInt blocksize)
34d71ae5a4SJacob Faibussowitsch {
35360ee056SFande Kong PetscInt subrstart, subrend, subrowsize, subcolsize, subcstart, subcend, rowsize, colsize;
368a2c336bSFande Kong PetscInt subrow, row, nz, *d_nnz, *o_nnz, i, j, dnz, onz, max_nz, *indices;
378a2c336bSFande Kong const PetscInt *idx;
388a2c336bSFande Kong const PetscScalar *values;
398a2c336bSFande Kong MPI_Comm comm;
40360ee056SFande Kong
41360ee056SFande Kong PetscFunctionBegin;
429566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)subinterp, &comm));
439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(subinterp, &subrstart, &subrend));
448a2c336bSFande Kong subrowsize = subrend - subrstart;
45360ee056SFande Kong rowsize = subrowsize * blocksize;
469566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(rowsize, &d_nnz, rowsize, &o_nnz));
479566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(subinterp, &subcstart, &subcend));
488a2c336bSFande Kong subcolsize = subcend - subcstart;
498a2c336bSFande Kong colsize = subcolsize * blocksize;
50360ee056SFande Kong max_nz = 0;
518a2c336bSFande Kong for (subrow = subrstart; subrow < subrend; subrow++) {
529566063dSJacob Faibussowitsch PetscCall(MatGetRow(subinterp, subrow, &nz, &idx, NULL));
53360ee056SFande Kong if (max_nz < nz) max_nz = nz;
549371c9d4SSatish Balay dnz = 0;
559371c9d4SSatish Balay onz = 0;
56360ee056SFande Kong for (i = 0; i < nz; i++) {
578a2c336bSFande Kong if (idx[i] >= subcstart && idx[i] < subcend) dnz++;
588a2c336bSFande Kong else onz++;
59360ee056SFande Kong }
60360ee056SFande Kong for (i = 0; i < blocksize; i++) {
61360ee056SFande Kong d_nnz[(subrow - subrstart) * blocksize + i] = dnz;
62360ee056SFande Kong o_nnz[(subrow - subrstart) * blocksize + i] = onz;
63360ee056SFande Kong }
649566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(subinterp, subrow, &nz, &idx, NULL));
65360ee056SFande Kong }
669566063dSJacob Faibussowitsch PetscCall(MatCreateAIJ(comm, rowsize, colsize, PETSC_DETERMINE, PETSC_DETERMINE, 0, d_nnz, 0, o_nnz, interp));
679566063dSJacob Faibussowitsch PetscCall(MatSetOption(*interp, MAT_IGNORE_OFF_PROC_ENTRIES, PETSC_TRUE));
689566063dSJacob Faibussowitsch PetscCall(MatSetOption(*interp, MAT_IGNORE_ZERO_ENTRIES, PETSC_TRUE));
699566063dSJacob Faibussowitsch PetscCall(MatSetOption(*interp, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE));
709566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(*interp));
71360ee056SFande Kong
729566063dSJacob Faibussowitsch PetscCall(MatSetUp(*interp));
739566063dSJacob Faibussowitsch PetscCall(PetscFree2(d_nnz, o_nnz));
749566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(max_nz, &indices));
758a2c336bSFande Kong for (subrow = subrstart; subrow < subrend; subrow++) {
769566063dSJacob Faibussowitsch PetscCall(MatGetRow(subinterp, subrow, &nz, &idx, &values));
77360ee056SFande Kong for (i = 0; i < blocksize; i++) {
78360ee056SFande Kong row = subrow * blocksize + i;
79ad540459SPierre Jolivet for (j = 0; j < nz; j++) indices[j] = idx[j] * blocksize + i;
809566063dSJacob Faibussowitsch PetscCall(MatSetValues(*interp, 1, &row, nz, indices, values, INSERT_VALUES));
81360ee056SFande Kong }
829566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(subinterp, subrow, &nz, &idx, &values));
83360ee056SFande Kong }
849566063dSJacob Faibussowitsch PetscCall(PetscFree(indices));
859566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*interp, MAT_FINAL_ASSEMBLY));
869566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*interp, MAT_FINAL_ASSEMBLY));
873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
88360ee056SFande Kong }
89360ee056SFande Kong
PCSetUp_HMG(PC pc)9066976f2fSJacob Faibussowitsch static PetscErrorCode PCSetUp_HMG(PC pc)
91d71ae5a4SJacob Faibussowitsch {
92360ee056SFande Kong Mat PA, submat;
93360ee056SFande Kong PC_MG *mg = (PC_MG *)pc->data;
94360ee056SFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx;
95360ee056SFande Kong MPI_Comm comm;
96360ee056SFande Kong PetscInt level;
97360ee056SFande Kong PetscInt num_levels;
988a2c336bSFande Kong Mat *operators, *interpolations;
998a2c336bSFande Kong PetscInt blocksize;
100fd2dd295SFande Kong const char *prefix;
10107a4832bSFande Kong PCMGGalerkinType galerkin;
102360ee056SFande Kong
103360ee056SFande Kong PetscFunctionBegin;
1049566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm));
105360ee056SFande Kong if (pc->setupcalled) {
1065f36d8f2SFande Kong if (hmg->reuseinterp) {
1075f36d8f2SFande Kong /* If we did not use Galerkin in the last call or we have a different sparsity pattern now,
1085f36d8f2SFande Kong * we have to build from scratch
10907a4832bSFande Kong * */
1109566063dSJacob Faibussowitsch PetscCall(PCMGGetGalerkin(pc, &galerkin));
1115f36d8f2SFande Kong if (galerkin == PC_MG_GALERKIN_NONE || pc->flag != SAME_NONZERO_PATTERN) pc->setupcalled = PETSC_FALSE;
1129566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, PC_MG_GALERKIN_PMAT));
1139566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc));
1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
115360ee056SFande Kong } else {
1169566063dSJacob Faibussowitsch PetscCall(PCReset_MG(pc));
117360ee056SFande Kong pc->setupcalled = PETSC_FALSE;
118360ee056SFande Kong }
119360ee056SFande Kong }
120360ee056SFande Kong
1218a2c336bSFande Kong /* Create an inner PC (GAMG or HYPRE) */
1228a2c336bSFande Kong if (!hmg->innerpc) {
1239566063dSJacob Faibussowitsch PetscCall(PCCreate(comm, &hmg->innerpc));
124fd2dd295SFande Kong /* If users do not set an inner pc type, we need to set a default value */
125fd2dd295SFande Kong if (!hmg->innerpctype) {
126fd2dd295SFande Kong /* If hypre is available, use hypre, otherwise, use gamg */
127a51cd166SPierre Jolivet #if PetscDefined(HAVE_HYPRE)
128f4f49eeaSPierre Jolivet PetscCall(PetscStrallocpy(PCHYPRE, &hmg->innerpctype));
129fd2dd295SFande Kong #else
130f4f49eeaSPierre Jolivet PetscCall(PetscStrallocpy(PCGAMG, &hmg->innerpctype));
131fd2dd295SFande Kong #endif
1328a2c336bSFande Kong }
1339566063dSJacob Faibussowitsch PetscCall(PCSetType(hmg->innerpc, hmg->innerpctype));
134360ee056SFande Kong }
1359566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &PA));
1368a2c336bSFande Kong /* Users need to correctly set a block size of matrix in order to use subspace coarsening */
1379566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(PA, &blocksize));
1388a2c336bSFande Kong if (blocksize <= 1) hmg->subcoarsening = PETSC_FALSE;
1398a2c336bSFande Kong /* Extract a submatrix for constructing subinterpolations */
1408a2c336bSFande Kong if (hmg->subcoarsening) {
1419566063dSJacob Faibussowitsch PetscCall(PCHMGExtractSubMatrix_Private(PA, &submat, MAT_INITIAL_MATRIX, hmg->component, blocksize));
142360ee056SFande Kong PA = submat;
143360ee056SFande Kong }
1449566063dSJacob Faibussowitsch PetscCall(PCSetOperators(hmg->innerpc, PA, PA));
14548a46eb9SPierre Jolivet if (hmg->subcoarsening) PetscCall(MatDestroy(&PA));
1468a2c336bSFande Kong /* Setup inner PC correctly. During this step, matrix will be coarsened */
1479566063dSJacob Faibussowitsch PetscCall(PCSetUseAmat(hmg->innerpc, PETSC_FALSE));
1489566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, &prefix));
1499566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)hmg->innerpc, prefix));
1509566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)hmg->innerpc, "hmg_inner_"));
1519566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(hmg->innerpc));
1529566063dSJacob Faibussowitsch PetscCall(PCSetUp(hmg->innerpc));
153360ee056SFande Kong
1548a2c336bSFande Kong /* Obtain interpolations IN PLACE. For BoomerAMG, (I,J,data) is reused to avoid memory overhead */
1559566063dSJacob Faibussowitsch PetscCall(PCGetInterpolations(hmg->innerpc, &num_levels, &interpolations));
1568a2c336bSFande Kong /* We can reuse the coarse operators when we do the full space coarsening */
15748a46eb9SPierre Jolivet if (!hmg->subcoarsening) PetscCall(PCGetCoarseOperators(hmg->innerpc, &num_levels, &operators));
158360ee056SFande Kong
1599566063dSJacob Faibussowitsch PetscCall(PCDestroy(&hmg->innerpc));
1600a545947SLisandro Dalcin hmg->innerpc = NULL;
1619566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels_MG(pc, num_levels, NULL));
1628a2c336bSFande Kong /* Set coarse matrices and interpolations to PCMG */
163360ee056SFande Kong for (level = num_levels - 1; level > 0; level--) {
1640a545947SLisandro Dalcin Mat P = NULL, pmat = NULL;
165360ee056SFande Kong Vec b, x, r;
1668a2c336bSFande Kong if (hmg->subcoarsening) {
1674bb91820SFande Kong if (hmg->usematmaij) {
1689566063dSJacob Faibussowitsch PetscCall(MatCreateMAIJ(interpolations[level - 1], blocksize, &P));
1699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&interpolations[level - 1]));
1704bb91820SFande Kong } else {
1718a2c336bSFande Kong /* Grow interpolation. In the future, we should use MAIJ */
1729566063dSJacob Faibussowitsch PetscCall(PCHMGExpandInterpolation_Private(interpolations[level - 1], &P, blocksize));
1739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&interpolations[level - 1]));
1744bb91820SFande Kong }
175360ee056SFande Kong } else {
1768a2c336bSFande Kong P = interpolations[level - 1];
177360ee056SFande Kong }
1789566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(P, &b, &r));
1799566063dSJacob Faibussowitsch PetscCall(PCMGSetInterpolation(pc, level, P));
1809566063dSJacob Faibussowitsch PetscCall(PCMGSetRestriction(pc, level, P));
1819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&P));
1828a2c336bSFande Kong /* We reuse the matrices when we do not do subspace coarsening */
1838a2c336bSFande Kong if ((level - 1) >= 0 && !hmg->subcoarsening) {
1848a2c336bSFande Kong pmat = operators[level - 1];
1859566063dSJacob Faibussowitsch PetscCall(PCMGSetOperators(pc, level - 1, pmat, pmat));
1869566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pmat));
187360ee056SFande Kong }
1889566063dSJacob Faibussowitsch PetscCall(PCMGSetRhs(pc, level - 1, b));
189360ee056SFande Kong
1909566063dSJacob Faibussowitsch PetscCall(PCMGSetR(pc, level, r));
1919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&r));
192360ee056SFande Kong
1939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b, &x));
1949566063dSJacob Faibussowitsch PetscCall(PCMGSetX(pc, level - 1, x));
1959566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
1969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
197360ee056SFande Kong }
1989566063dSJacob Faibussowitsch PetscCall(PetscFree(interpolations));
19948a46eb9SPierre Jolivet if (!hmg->subcoarsening) PetscCall(PetscFree(operators));
2008a2c336bSFande Kong /* Turn Galerkin off when we already have coarse operators */
2019566063dSJacob Faibussowitsch PetscCall(PCMGSetGalerkin(pc, hmg->subcoarsening ? PC_MG_GALERKIN_PMAT : PC_MG_GALERKIN_NONE));
2029566063dSJacob Faibussowitsch PetscCall(PCSetDM(pc, NULL));
2039566063dSJacob Faibussowitsch PetscCall(PCSetUseAmat(pc, PETSC_FALSE));
204d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)pc);
205dbbe0bcdSBarry Smith PetscCall(PCSetFromOptions_MG(pc, PetscOptionsObject)); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */
206d0609cedSBarry Smith PetscOptionsEnd();
2079566063dSJacob Faibussowitsch PetscCall(PCSetUp_MG(pc));
2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
209360ee056SFande Kong }
210360ee056SFande Kong
PCDestroy_HMG(PC pc)21166976f2fSJacob Faibussowitsch static PetscErrorCode PCDestroy_HMG(PC pc)
212d71ae5a4SJacob Faibussowitsch {
213360ee056SFande Kong PC_MG *mg = (PC_MG *)pc->data;
2148a2c336bSFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx;
215360ee056SFande Kong
216360ee056SFande Kong PetscFunctionBegin;
2179566063dSJacob Faibussowitsch PetscCall(PCDestroy(&hmg->innerpc));
2189566063dSJacob Faibussowitsch PetscCall(PetscFree(hmg->innerpctype));
2199566063dSJacob Faibussowitsch PetscCall(PetscFree(hmg));
2209566063dSJacob Faibussowitsch PetscCall(PCDestroy_MG(pc));
221fd2dd295SFande Kong
2229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetReuseInterpolation_C", NULL));
2239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetUseSubspaceCoarsening_C", NULL));
2249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetInnerPCType_C", NULL));
2259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetCoarseningComponent_C", NULL));
2262e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGUseMatMAIJ_C", NULL));
2273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
228360ee056SFande Kong }
229360ee056SFande Kong
PCView_HMG(PC pc,PetscViewer viewer)23066976f2fSJacob Faibussowitsch static PetscErrorCode PCView_HMG(PC pc, PetscViewer viewer)
231d71ae5a4SJacob Faibussowitsch {
232360ee056SFande Kong PC_MG *mg = (PC_MG *)pc->data;
2338a2c336bSFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx;
234*9f196a02SMartin Diehl PetscBool isascii;
235360ee056SFande Kong
236360ee056SFande Kong PetscFunctionBegin;
237*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
238*9f196a02SMartin Diehl if (isascii) {
2399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Reuse interpolation: %s\n", hmg->reuseinterp ? "true" : "false"));
2409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Use subspace coarsening: %s\n", hmg->subcoarsening ? "true" : "false"));
24163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Coarsening component: %" PetscInt_FMT " \n", hmg->component));
2429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Use MatMAIJ: %s \n", hmg->usematmaij ? "true" : "false"));
2439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Inner PC type: %s \n", hmg->innerpctype));
244360ee056SFande Kong }
2459566063dSJacob Faibussowitsch PetscCall(PCView_MG(pc, viewer));
2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
247360ee056SFande Kong }
248360ee056SFande Kong
PCSetFromOptions_HMG(PC pc,PetscOptionItems PetscOptionsObject)249ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_HMG(PC pc, PetscOptionItems PetscOptionsObject)
250d71ae5a4SJacob Faibussowitsch {
251360ee056SFande Kong PC_MG *mg = (PC_MG *)pc->data;
2528a2c336bSFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx;
253360ee056SFande Kong
254360ee056SFande Kong PetscFunctionBegin;
255d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "HMG");
2569566063dSJacob 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));
2579566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_hmg_use_subspace_coarsening", "Use the subspace coarsening to compute the interpolations", "PCHMGSetUseSubspaceCoarsening", hmg->subcoarsening, &hmg->subcoarsening, NULL));
2589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_hmg_use_matmaij", "Use MatMAIJ store interpolation for saving memory", "PCHMGSetInnerPCType", hmg->usematmaij, &hmg->usematmaij, NULL));
2599566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_hmg_coarsening_component", "Which component is chosen for the subspace-based coarsening algorithm", "PCHMGSetCoarseningComponent", hmg->component, &hmg->component, NULL));
260d0609cedSBarry Smith PetscOptionsHeadEnd();
2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
262360ee056SFande Kong }
263360ee056SFande Kong
PCHMGSetReuseInterpolation_HMG(PC pc,PetscBool reuse)264d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGSetReuseInterpolation_HMG(PC pc, PetscBool reuse)
265d71ae5a4SJacob Faibussowitsch {
26607a4832bSFande Kong PC_MG *mg = (PC_MG *)pc->data;
26707a4832bSFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx;
268fd2dd295SFande Kong
269fd2dd295SFande Kong PetscFunctionBegin;
270fd2dd295SFande Kong hmg->reuseinterp = reuse;
2713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
272fd2dd295SFande Kong }
273fd2dd295SFande Kong
274b155ba7fSFande Kong /*@
275f1580f4eSBarry Smith PCHMGSetReuseInterpolation - Reuse the interpolation matrices in `PCHMG` after changing the matrices numerical values
2768a2c336bSFande Kong
277c3339decSBarry Smith Logically Collective
2788a2c336bSFande Kong
2798a2c336bSFande Kong Input Parameters:
280f1580f4eSBarry Smith + pc - the `PCHMG` context
281f1580f4eSBarry Smith - reuse - `PETSC_TRUE` indicates that `PCHMG` will reuse the interpolations
2828a2c336bSFande Kong
283f1580f4eSBarry Smith Options Database Key:
28449c604d5SFande Kong . -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
2858a2c336bSFande Kong
2868a2c336bSFande Kong Level: beginner
2878a2c336bSFande Kong
2880b4b7b1cSBarry Smith Note:
2890b4b7b1cSBarry Smith 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`
2900b4b7b1cSBarry Smith
291562efe2eSBarry Smith .seealso: [](ch_ksp), `PCHMG`, `PCGAMG`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetCoarseningComponent()`, `PCHMGSetInnerPCType()`
292b155ba7fSFande Kong @*/
PCHMGSetReuseInterpolation(PC pc,PetscBool reuse)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
PCHMGSetUseSubspaceCoarsening_HMG(PC pc,PetscBool subspace)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 @*/
PCHMGSetUseSubspaceCoarsening(PC pc,PetscBool subspace)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
PCHMGSetInnerPCType_HMG(PC pc,PCType type)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;
341f4f49eeaSPierre Jolivet PetscCall(PetscStrallocpy(type, &hmg->innerpctype));
3423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3438a2c336bSFande Kong }
3448a2c336bSFande Kong
3455d83a8b1SBarry Smith /*@
3460b4b7b1cSBarry Smith PCHMGSetInnerPCType - Set an inner `PC` type to be used in the `PCHMG` preconditioner. That is the method used to compute
3470b4b7b1cSBarry Smith the hierarchy of restriction operators.
3488a2c336bSFande Kong
349c3339decSBarry Smith Logically Collective
3508a2c336bSFande Kong
3518a2c336bSFande Kong Input Parameters:
352f1580f4eSBarry Smith + pc - the `PCHMG` context
353f1580f4eSBarry Smith - type - `PCHYPRE` or `PCGAMG` coarsening algorithm
3548a2c336bSFande Kong
355f1580f4eSBarry Smith Options Database Key:
35649c604d5SFande Kong . -hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix
3578a2c336bSFande Kong
3588a2c336bSFande Kong Level: beginner
3598a2c336bSFande Kong
360562efe2eSBarry Smith .seealso: [](ch_ksp), `PCHMG`, `PCType`, `PCHMGSetReuseInterpolation()`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetCoarseningComponent()`
361b155ba7fSFande Kong @*/
PCHMGSetInnerPCType(PC pc,PCType type)362d71ae5a4SJacob Faibussowitsch PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type)
363d71ae5a4SJacob Faibussowitsch {
3648a2c336bSFande Kong PetscFunctionBegin;
3658a2c336bSFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
366cac4c232SBarry Smith PetscUseMethod(pc, "PCHMGSetInnerPCType_C", (PC, PCType), (pc, type));
3673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
3688a2c336bSFande Kong }
3698a2c336bSFande Kong
PCHMGSetCoarseningComponent_HMG(PC pc,PetscInt component)370d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGSetCoarseningComponent_HMG(PC pc, PetscInt component)
371d71ae5a4SJacob Faibussowitsch {
37249c604d5SFande Kong PC_MG *mg = (PC_MG *)pc->data;
37349c604d5SFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx;
37449c604d5SFande Kong
37549c604d5SFande Kong PetscFunctionBegin;
37649c604d5SFande Kong hmg->component = component;
3773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
37849c604d5SFande Kong }
37949c604d5SFande Kong
380b155ba7fSFande Kong /*@
3810b4b7b1cSBarry Smith PCHMGSetCoarseningComponent - Set which component of the PDE is used for the subspace-based coarsening algorithm in the preconditioner `PCHMG`
38249c604d5SFande Kong
383c3339decSBarry Smith Logically Collective
38449c604d5SFande Kong
38549c604d5SFande Kong Input Parameters:
386f1580f4eSBarry Smith + pc - the `PCHMG` context
387f1580f4eSBarry Smith - component - which component `PC` will coarsen
38849c604d5SFande Kong
389f1580f4eSBarry Smith Options Database Key:
390f1580f4eSBarry Smith . -pc_hmg_coarsening_component <i> - Which component is chosen for the subspace-based coarsening algorithm
39149c604d5SFande Kong
39249c604d5SFande Kong Level: beginner
39349c604d5SFande Kong
3940b4b7b1cSBarry Smith Note:
3950b4b7b1cSBarry Smith By default it uses the first component
3960b4b7b1cSBarry Smith
397562efe2eSBarry Smith .seealso: [](ch_ksp), `PCHMG`, `PCType`, `PCGAMG`, `PCHMGSetReuseInterpolation()`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetInnerPCType()`
398b155ba7fSFande Kong @*/
PCHMGSetCoarseningComponent(PC pc,PetscInt component)399d71ae5a4SJacob Faibussowitsch PetscErrorCode PCHMGSetCoarseningComponent(PC pc, PetscInt component)
400d71ae5a4SJacob Faibussowitsch {
40149c604d5SFande Kong PetscFunctionBegin;
40249c604d5SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
403cac4c232SBarry Smith PetscUseMethod(pc, "PCHMGSetCoarseningComponent_C", (PC, PetscInt), (pc, component));
4043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
40549c604d5SFande Kong }
40649c604d5SFande Kong
PCHMGUseMatMAIJ_HMG(PC pc,PetscBool usematmaij)407d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCHMGUseMatMAIJ_HMG(PC pc, PetscBool usematmaij)
408d71ae5a4SJacob Faibussowitsch {
40949c604d5SFande Kong PC_MG *mg = (PC_MG *)pc->data;
41049c604d5SFande Kong PC_HMG *hmg = (PC_HMG *)mg->innerctx;
41149c604d5SFande Kong
41249c604d5SFande Kong PetscFunctionBegin;
41349c604d5SFande Kong hmg->usematmaij = usematmaij;
4143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
41549c604d5SFande Kong }
41649c604d5SFande Kong
417b155ba7fSFande Kong /*@
4180b4b7b1cSBarry Smith PCHMGUseMatMAIJ - Set a flag that indicates if or not to use `MATMAIJ` for the interpolation matrices to save memory
41949c604d5SFande Kong
420c3339decSBarry Smith Logically Collective
42149c604d5SFande Kong
42249c604d5SFande Kong Input Parameters:
423f1580f4eSBarry Smith + pc - the `PCHMG` context
424f1580f4eSBarry Smith - usematmaij - `PETSC_TRUE` (default) to use `MATMAIJ` for interpolations.
42549c604d5SFande Kong
426f1580f4eSBarry Smith Options Database Key:
42749c604d5SFande Kong . -pc_hmg_use_matmaij - <true | false >
42849c604d5SFande Kong
42949c604d5SFande Kong Level: beginner
43049c604d5SFande Kong
431562efe2eSBarry Smith .seealso: [](ch_ksp), `PCHMG`, `PCType`, `PCGAMG`
432b155ba7fSFande Kong @*/
PCHMGUseMatMAIJ(PC pc,PetscBool usematmaij)433d71ae5a4SJacob Faibussowitsch PetscErrorCode PCHMGUseMatMAIJ(PC pc, PetscBool usematmaij)
434d71ae5a4SJacob Faibussowitsch {
43549c604d5SFande Kong PetscFunctionBegin;
43649c604d5SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
437cac4c232SBarry Smith PetscUseMethod(pc, "PCHMGUseMatMAIJ_C", (PC, PetscBool), (pc, usematmaij));
4383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
43949c604d5SFande Kong }
44049c604d5SFande Kong
4418a2c336bSFande Kong /*MC
4420b4b7b1cSBarry Smith PCHMG - Preconditioner for multiple component PDE problems that constructs a hierarchy of restriction operators to coarse grid problems using the submatrix of
4430b4b7b1cSBarry Smith 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`
4440b4b7b1cSBarry Smith multigrid preconditioner. This results in a much more efficient to build and apply preconditioner than using `PCGAMG` on the entire system {cite}`kong2020highly`.
445360ee056SFande Kong
446360ee056SFande Kong Options Database Keys:
4470b4b7b1cSBarry Smith + -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.
4480b4b7b1cSBarry Smith . -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).
4490b4b7b1cSBarry Smith . -hmg_inner_pc_type <hypre, gamg, ...> - What method to use to generate the hierarchy of restriction operators
450f1580f4eSBarry Smith - -pc_hmg_use_matmaij <true | false> - Whether or not to use `MATMAIJ` for multicomponent problems for saving memory
451360ee056SFande Kong
452f1580f4eSBarry Smith Level: intermediate
453360ee056SFande Kong
454f1580f4eSBarry Smith Note:
455f1580f4eSBarry Smith `MatSetBlockSize()` must be called on the linear system matrix to set the number of components of the PDE.
456360ee056SFande Kong
457562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMG`, `PCHYPRE`, `PCHMG`, `PCGetCoarseOperators()`, `PCGetInterpolations()`,
4580b4b7b1cSBarry Smith `PCHMGSetReuseInterpolation()`, `PCHMGSetUseSubspaceCoarsening()`, `PCHMGSetInnerPCType()`, `PCGAMG`
459360ee056SFande Kong M*/
PCCreate_HMG(PC pc)460d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc)
461d71ae5a4SJacob Faibussowitsch {
462360ee056SFande Kong PC_HMG *hmg;
463360ee056SFande Kong PC_MG *mg;
464360ee056SFande Kong
465360ee056SFande Kong PetscFunctionBegin;
466360ee056SFande Kong /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
467dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, destroy);
4680a545947SLisandro Dalcin pc->data = NULL;
4699566063dSJacob Faibussowitsch PetscCall(PetscFree(((PetscObject)pc)->type_name));
470360ee056SFande Kong
4719566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCMG));
4729566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCHMG));
4739566063dSJacob Faibussowitsch PetscCall(PetscNew(&hmg));
474360ee056SFande Kong
475360ee056SFande Kong mg = (PC_MG *)pc->data;
476360ee056SFande Kong mg->innerctx = hmg;
477360ee056SFande Kong hmg->reuseinterp = PETSC_FALSE;
4788a2c336bSFande Kong hmg->subcoarsening = PETSC_FALSE;
4794bb91820SFande Kong hmg->usematmaij = PETSC_TRUE;
48049c604d5SFande Kong hmg->component = 0;
4818a2c336bSFande Kong hmg->innerpc = NULL;
482360ee056SFande Kong
483360ee056SFande Kong pc->ops->setfromoptions = PCSetFromOptions_HMG;
484360ee056SFande Kong pc->ops->view = PCView_HMG;
485360ee056SFande Kong pc->ops->destroy = PCDestroy_HMG;
486360ee056SFande Kong pc->ops->setup = PCSetUp_HMG;
487fd2dd295SFande Kong
4889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetReuseInterpolation_C", PCHMGSetReuseInterpolation_HMG));
4899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetUseSubspaceCoarsening_C", PCHMGSetUseSubspaceCoarsening_HMG));
4909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetInnerPCType_C", PCHMGSetInnerPCType_HMG));
4919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGSetCoarseningComponent_C", PCHMGSetCoarseningComponent_HMG));
4929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCHMGUseMatMAIJ_C", PCHMGUseMatMAIJ_HMG));
4933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
494360ee056SFande Kong }
495