xref: /petsc/src/ksp/pc/impls/hmg/hmg.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
1360ee056SFande Kong #include <petscdm.h>
2360ee056SFande Kong #include <petscctable.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 
16360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems*,PC);
17360ee056SFande Kong PetscErrorCode PCReset_MG(PC);
188a2c336bSFande Kong 
1949c604d5SFande Kong static PetscErrorCode PCHMGExtractSubMatrix_Private(Mat pmat,Mat *submat,MatReuse reuse,PetscInt component,PetscInt blocksize)
20360ee056SFande Kong {
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));
272c71b3e2SJacob Faibussowitsch   PetscCheckFalse(component>=blocksize,comm,PETSC_ERR_ARG_INCOMP,"Component %D should be less than block size %D ",component,blocksize);
289566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(pmat,&rstart,&rend));
292c71b3e2SJacob Faibussowitsch   PetscCheckFalse((rend-rstart)%blocksize != 0,comm,PETSC_ERR_ARG_INCOMP,"Block size %D is inconsistent for [%D, %D) ",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));
33360ee056SFande Kong   PetscFunctionReturn(0);
34360ee056SFande Kong }
35360ee056SFande Kong 
368a2c336bSFande Kong static PetscErrorCode PCHMGExpandInterpolation_Private(Mat subinterp, Mat *interp, PetscInt blocksize)
37360ee056SFande Kong {
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;
57360ee056SFande Kong     dnz = 0; onz = 0;
58360ee056SFande Kong     for (i=0;i<nz;i++) {
598a2c336bSFande Kong       if (idx[i]>=subcstart && idx[i]<subcend) dnz++;
608a2c336bSFande Kong       else onz++;
61360ee056SFande Kong     }
62360ee056SFande Kong     for (i=0;i<blocksize;i++) {
63360ee056SFande Kong       d_nnz[(subrow-subrstart)*blocksize+i] = dnz;
64360ee056SFande Kong       o_nnz[(subrow-subrstart)*blocksize+i] = onz;
65360ee056SFande Kong     }
669566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(subinterp,subrow,&nz,&idx,NULL));
67360ee056SFande Kong   }
689566063dSJacob Faibussowitsch   PetscCall(MatCreateAIJ(comm,rowsize,colsize,PETSC_DETERMINE,PETSC_DETERMINE,0,d_nnz,0,o_nnz,interp));
699566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*interp,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE));
709566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*interp,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE));
719566063dSJacob Faibussowitsch   PetscCall(MatSetOption(*interp,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
729566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(*interp));
73360ee056SFande Kong 
749566063dSJacob Faibussowitsch   PetscCall(MatSetUp(*interp));
759566063dSJacob Faibussowitsch   PetscCall(PetscFree2(d_nnz,o_nnz));
769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(max_nz,&indices));
778a2c336bSFande Kong   for (subrow=subrstart; subrow<subrend; subrow++) {
789566063dSJacob Faibussowitsch     PetscCall(MatGetRow(subinterp,subrow,&nz,&idx,&values));
79360ee056SFande Kong     for (i=0;i<blocksize;i++) {
80360ee056SFande Kong       row = subrow*blocksize+i;
81360ee056SFande Kong       for (j=0;j<nz;j++) {
82360ee056SFande Kong         indices[j] = idx[j]*blocksize+i;
83360ee056SFande Kong       }
849566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*interp,1,&row,nz,indices,values,INSERT_VALUES));
85360ee056SFande Kong     }
869566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(subinterp,subrow,&nz,&idx,&values));
87360ee056SFande Kong   }
889566063dSJacob Faibussowitsch   PetscCall(PetscFree(indices));
899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*interp,MAT_FINAL_ASSEMBLY));
909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*interp,MAT_FINAL_ASSEMBLY));
91360ee056SFande Kong   PetscFunctionReturn(0);
92360ee056SFande Kong }
93360ee056SFande Kong 
94360ee056SFande Kong PetscErrorCode PCSetUp_HMG(PC pc)
95360ee056SFande Kong {
96360ee056SFande Kong   PetscErrorCode     ierr;
97360ee056SFande Kong   Mat                PA, submat;
98360ee056SFande Kong   PC_MG              *mg   = (PC_MG*)pc->data;
99360ee056SFande Kong   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
100360ee056SFande Kong   MPI_Comm           comm;
101360ee056SFande Kong   PetscInt           level;
102360ee056SFande Kong   PetscInt           num_levels;
1038a2c336bSFande Kong   Mat                *operators,*interpolations;
1048a2c336bSFande Kong   PetscInt           blocksize;
105fd2dd295SFande Kong   const char         *prefix;
10607a4832bSFande Kong   PCMGGalerkinType   galerkin;
107360ee056SFande Kong 
108360ee056SFande Kong   PetscFunctionBegin;
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)pc,&comm));
110360ee056SFande Kong   if (pc->setupcalled) {
1115f36d8f2SFande Kong    if (hmg->reuseinterp) {
1125f36d8f2SFande Kong      /* If we did not use Galerkin in the last call or we have a different sparsity pattern now,
1135f36d8f2SFande Kong       * we have to build from scratch
11407a4832bSFande Kong       * */
1159566063dSJacob Faibussowitsch      PetscCall(PCMGGetGalerkin(pc,&galerkin));
1165f36d8f2SFande Kong      if (galerkin == PC_MG_GALERKIN_NONE || pc->flag != SAME_NONZERO_PATTERN) pc->setupcalled = PETSC_FALSE;
1179566063dSJacob Faibussowitsch      PetscCall(PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT));
1189566063dSJacob Faibussowitsch      PetscCall(PCSetUp_MG(pc));
119360ee056SFande Kong      PetscFunctionReturn(0);
120360ee056SFande Kong     } else {
1219566063dSJacob Faibussowitsch      PetscCall(PCReset_MG(pc));
122360ee056SFande Kong      pc->setupcalled = PETSC_FALSE;
123360ee056SFande Kong     }
124360ee056SFande Kong   }
125360ee056SFande Kong 
1268a2c336bSFande Kong   /* Create an inner PC (GAMG or HYPRE) */
1278a2c336bSFande Kong   if (!hmg->innerpc) {
1289566063dSJacob Faibussowitsch     PetscCall(PCCreate(comm,&hmg->innerpc));
129fd2dd295SFande Kong     /* If users do not set an inner pc type, we need to set a default value */
130fd2dd295SFande Kong     if (!hmg->innerpctype) {
131fd2dd295SFande Kong     /* If hypre is available, use hypre, otherwise, use gamg */
132fd2dd295SFande Kong #if PETSC_HAVE_HYPRE
1339566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(PCHYPRE,&(hmg->innerpctype)));
134fd2dd295SFande Kong #else
1359566063dSJacob Faibussowitsch       PetscCall(PetscStrallocpy(PCGAMG,&(hmg->innerpctype)));
136fd2dd295SFande Kong #endif
1378a2c336bSFande Kong     }
1389566063dSJacob Faibussowitsch     PetscCall(PCSetType(hmg->innerpc,hmg->innerpctype));
139360ee056SFande Kong   }
1409566063dSJacob Faibussowitsch   PetscCall(PCGetOperators(pc,NULL,&PA));
1418a2c336bSFande Kong   /* Users need to correctly set a block size of matrix in order to use subspace coarsening */
1429566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(PA,&blocksize));
1438a2c336bSFande Kong   if (blocksize<=1) hmg->subcoarsening = PETSC_FALSE;
1448a2c336bSFande Kong   /* Extract a submatrix for constructing subinterpolations */
1458a2c336bSFande Kong   if (hmg->subcoarsening) {
1469566063dSJacob Faibussowitsch     PetscCall(PCHMGExtractSubMatrix_Private(PA,&submat,MAT_INITIAL_MATRIX,hmg->component,blocksize));
147360ee056SFande Kong     PA = submat;
148360ee056SFande Kong   }
1499566063dSJacob Faibussowitsch   PetscCall(PCSetOperators(hmg->innerpc,PA,PA));
1508a2c336bSFande Kong   if (hmg->subcoarsening) {
1519566063dSJacob Faibussowitsch    PetscCall(MatDestroy(&PA));
152360ee056SFande Kong   }
1538a2c336bSFande Kong   /* Setup inner PC correctly. During this step, matrix will be coarsened */
1549566063dSJacob Faibussowitsch   PetscCall(PCSetUseAmat(hmg->innerpc,PETSC_FALSE));
1559566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc,&prefix));
1569566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)hmg->innerpc,prefix));
1579566063dSJacob Faibussowitsch   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)hmg->innerpc,"hmg_inner_"));
1589566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions(hmg->innerpc));
1599566063dSJacob Faibussowitsch   PetscCall(PCSetUp(hmg->innerpc));
160360ee056SFande Kong 
1618a2c336bSFande Kong   /* Obtain interpolations IN PLACE. For BoomerAMG, (I,J,data) is reused to avoid memory overhead */
1629566063dSJacob Faibussowitsch   PetscCall(PCGetInterpolations(hmg->innerpc,&num_levels,&interpolations));
1638a2c336bSFande Kong   /* We can reuse the coarse operators when we do the full space coarsening */
1648a2c336bSFande Kong   if (!hmg->subcoarsening) {
1659566063dSJacob Faibussowitsch     PetscCall(PCGetCoarseOperators(hmg->innerpc,&num_levels,&operators));
166360ee056SFande Kong   }
167360ee056SFande Kong 
1689566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&hmg->innerpc));
1690a545947SLisandro Dalcin   hmg->innerpc = NULL;
1709566063dSJacob Faibussowitsch   PetscCall(PCMGSetLevels_MG(pc,num_levels,NULL));
1718a2c336bSFande Kong   /* Set coarse matrices and interpolations to PCMG */
172360ee056SFande Kong   for (level=num_levels-1; level>0; level--) {
1730a545947SLisandro Dalcin     Mat P=NULL, pmat=NULL;
174360ee056SFande Kong     Vec b, x,r;
1758a2c336bSFande Kong     if (hmg->subcoarsening) {
1764bb91820SFande Kong       if (hmg->usematmaij) {
1779566063dSJacob Faibussowitsch         PetscCall(MatCreateMAIJ(interpolations[level-1],blocksize,&P));
1789566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&interpolations[level-1]));
1794bb91820SFande Kong       } else {
1808a2c336bSFande Kong         /* Grow interpolation. In the future, we should use MAIJ */
1819566063dSJacob Faibussowitsch         PetscCall(PCHMGExpandInterpolation_Private(interpolations[level-1],&P,blocksize));
1829566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&interpolations[level-1]));
1834bb91820SFande Kong       }
184360ee056SFande Kong     } else {
1858a2c336bSFande Kong       P = interpolations[level-1];
186360ee056SFande Kong     }
1879566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(P,&b,&r));
1889566063dSJacob Faibussowitsch     PetscCall(PCMGSetInterpolation(pc,level,P));
1899566063dSJacob Faibussowitsch     PetscCall(PCMGSetRestriction(pc,level,P));
1909566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&P));
1918a2c336bSFande Kong     /* We reuse the matrices when we do not do subspace coarsening */
1928a2c336bSFande Kong     if ((level-1)>=0 && !hmg->subcoarsening) {
1938a2c336bSFande Kong       pmat = operators[level-1];
1949566063dSJacob Faibussowitsch       PetscCall(PCMGSetOperators(pc,level-1,pmat,pmat));
1959566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&pmat));
196360ee056SFande Kong     }
1979566063dSJacob Faibussowitsch     PetscCall(PCMGSetRhs(pc,level-1,b));
198360ee056SFande Kong 
1999566063dSJacob Faibussowitsch     PetscCall(PCMGSetR(pc,level,r));
2009566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&r));
201360ee056SFande Kong 
2029566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(b,&x));
2039566063dSJacob Faibussowitsch     PetscCall(PCMGSetX(pc,level-1,x));
2049566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&x));
2059566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&b));
206360ee056SFande Kong   }
2079566063dSJacob Faibussowitsch   PetscCall(PetscFree(interpolations));
2088a2c336bSFande Kong   if (!hmg->subcoarsening) {
2099566063dSJacob Faibussowitsch     PetscCall(PetscFree(operators));
2108a2c336bSFande Kong   }
2118a2c336bSFande Kong   /* Turn Galerkin off when we already have coarse operators */
2129566063dSJacob Faibussowitsch   PetscCall(PCMGSetGalerkin(pc,hmg->subcoarsening ? PC_MG_GALERKIN_PMAT:PC_MG_GALERKIN_NONE));
2139566063dSJacob Faibussowitsch   PetscCall(PCSetDM(pc,NULL));
2149566063dSJacob Faibussowitsch   PetscCall(PCSetUseAmat(pc,PETSC_FALSE));
2159566063dSJacob Faibussowitsch   ierr = PetscObjectOptionsBegin((PetscObject)pc);PetscCall(ierr);
2169566063dSJacob Faibussowitsch   PetscCall(PCSetFromOptions_MG(PetscOptionsObject,pc)); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */
2179566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
2189566063dSJacob Faibussowitsch   PetscCall(PCSetUp_MG(pc));
219360ee056SFande Kong   PetscFunctionReturn(0);
220360ee056SFande Kong }
221360ee056SFande Kong 
222360ee056SFande Kong PetscErrorCode PCDestroy_HMG(PC pc)
223360ee056SFande Kong {
224360ee056SFande Kong   PC_MG          *mg  = (PC_MG*)pc->data;
2258a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
226360ee056SFande Kong 
227360ee056SFande Kong   PetscFunctionBegin;
2289566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&hmg->innerpc));
2299566063dSJacob Faibussowitsch   PetscCall(PetscFree(hmg->innerpctype));
2309566063dSJacob Faibussowitsch   PetscCall(PetscFree(hmg));
2319566063dSJacob Faibussowitsch   PetscCall(PCDestroy_MG(pc));
232fd2dd295SFande Kong 
2339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",NULL));
2349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",NULL));
2359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",NULL));
2369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetCoarseningComponent_C",NULL));
237360ee056SFande Kong   PetscFunctionReturn(0);
238360ee056SFande Kong }
239360ee056SFande Kong 
240360ee056SFande Kong PetscErrorCode PCView_HMG(PC pc,PetscViewer viewer)
241360ee056SFande Kong {
242360ee056SFande Kong   PC_MG          *mg = (PC_MG*)pc->data;
2438a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
244360ee056SFande Kong   PetscBool      iascii;
245360ee056SFande Kong 
246360ee056SFande Kong   PetscFunctionBegin;
2479566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
248360ee056SFande Kong   if (iascii) {
2499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer," Reuse interpolation: %s\n",hmg->reuseinterp? "true":"false"));
2509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer," Use subspace coarsening: %s\n",hmg->subcoarsening? "true":"false"));
2519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer," Coarsening component: %D \n",hmg->component));
2529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer," Use MatMAIJ: %s \n",hmg->usematmaij? "true":"false"));
2539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer," Inner PC type: %s \n",hmg->innerpctype));
254360ee056SFande Kong   }
2559566063dSJacob Faibussowitsch   PetscCall(PCView_MG(pc,viewer));
256360ee056SFande Kong   PetscFunctionReturn(0);
257360ee056SFande Kong }
258360ee056SFande Kong 
259360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems *PetscOptionsObject,PC pc)
260360ee056SFande Kong {
261360ee056SFande Kong   PC_MG          *mg = (PC_MG*)pc->data;
2628a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
263360ee056SFande Kong 
264360ee056SFande Kong   PetscFunctionBegin;
2659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject,"HMG"));
2669566063dSJacob 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));
2679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_hmg_use_subspace_coarsening","Use the subspace coarsening to compute the interpolations","PCHMGSetUseSubspaceCoarsening",hmg->subcoarsening,&hmg->subcoarsening,NULL));
2689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_hmg_use_matmaij","Use MatMAIJ store interpolation for saving memory","PCHMGSetInnerPCType",hmg->usematmaij,&hmg->usematmaij,NULL));
2699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_hmg_coarsening_component","Which component is chosen for the subspace-based coarsening algorithm","PCHMGSetCoarseningComponent",hmg->component,&hmg->component,NULL));
2709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
271360ee056SFande Kong   PetscFunctionReturn(0);
272360ee056SFande Kong }
273360ee056SFande Kong 
274fd2dd295SFande Kong static PetscErrorCode PCHMGSetReuseInterpolation_HMG(PC pc, PetscBool reuse)
275fd2dd295SFande Kong {
27607a4832bSFande Kong   PC_MG          *mg  = (PC_MG*)pc->data;
27707a4832bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
278fd2dd295SFande Kong 
279fd2dd295SFande Kong   PetscFunctionBegin;
280fd2dd295SFande Kong   hmg->reuseinterp = reuse;
281fd2dd295SFande Kong   PetscFunctionReturn(0);
282fd2dd295SFande Kong }
283fd2dd295SFande Kong 
284b155ba7fSFande Kong /*@
2858a2c336bSFande Kong    PCHMGSetReuseInterpolation - Reuse interpolation matrices in HMG
2868a2c336bSFande Kong 
2878a2c336bSFande Kong    Logically Collective on PC
2888a2c336bSFande Kong 
2898a2c336bSFande Kong    Input Parameters:
2908a2c336bSFande Kong +  pc - the HMG context
2918a2c336bSFande Kong -  reuse - True indicates that HMG will reuse the interpolations
2928a2c336bSFande Kong 
2938a2c336bSFande Kong    Options Database Keys:
29449c604d5SFande Kong .  -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
2958a2c336bSFande Kong 
2968a2c336bSFande Kong    Level: beginner
2978a2c336bSFande Kong 
2988a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, reuse, set
2998a2c336bSFande Kong 
3008a2c336bSFande Kong .seealso: PCHMG
301b155ba7fSFande Kong @*/
3028a2c336bSFande Kong PetscErrorCode PCHMGSetReuseInterpolation(PC pc, PetscBool reuse)
3038a2c336bSFande Kong {
304fd2dd295SFande Kong   PetscFunctionBegin;
305fd2dd295SFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
306*cac4c232SBarry Smith   PetscUseMethod(pc,"PCHMGSetReuseInterpolation_C",(PC,PetscBool),(pc,reuse));
307fd2dd295SFande Kong   PetscFunctionReturn(0);
308fd2dd295SFande Kong }
309fd2dd295SFande Kong 
310fd2dd295SFande Kong static PetscErrorCode PCHMGSetUseSubspaceCoarsening_HMG(PC pc, PetscBool subspace)
311fd2dd295SFande Kong {
31207a4832bSFande Kong   PC_MG          *mg  = (PC_MG*)pc->data;
31307a4832bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
3148a2c336bSFande Kong 
3158a2c336bSFande Kong   PetscFunctionBegin;
316fd2dd295SFande Kong   hmg->subcoarsening = subspace;
3178a2c336bSFande Kong   PetscFunctionReturn(0);
3188a2c336bSFande Kong }
3198a2c336bSFande Kong 
320b155ba7fSFande Kong /*@
3218a2c336bSFande Kong    PCHMGSetUseSubspaceCoarsening - Use subspace coarsening in HMG
3228a2c336bSFande Kong 
3238a2c336bSFande Kong    Logically Collective on PC
3248a2c336bSFande Kong 
3258a2c336bSFande Kong    Input Parameters:
3268a2c336bSFande Kong +  pc - the HMG context
3278a2c336bSFande Kong -  reuse - True indicates that HMG will use the subspace coarsening
3288a2c336bSFande Kong 
3298a2c336bSFande Kong    Options Database Keys:
33049c604d5SFande Kong .  -pc_hmg_use_subspace_coarsening  <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix).
3318a2c336bSFande Kong 
3328a2c336bSFande Kong    Level: beginner
3338a2c336bSFande Kong 
3348a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, subspace, coarsening
3358a2c336bSFande Kong 
3368a2c336bSFande Kong .seealso: PCHMG
337b155ba7fSFande Kong @*/
3388a2c336bSFande Kong PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC pc, PetscBool subspace)
3398a2c336bSFande Kong {
3408a2c336bSFande Kong   PetscFunctionBegin;
3418a2c336bSFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
342*cac4c232SBarry Smith   PetscUseMethod(pc,"PCHMGSetUseSubspaceCoarsening_C",(PC,PetscBool),(pc,subspace));
343fd2dd295SFande Kong   PetscFunctionReturn(0);
344fd2dd295SFande Kong }
345fd2dd295SFande Kong 
346fd2dd295SFande Kong static PetscErrorCode PCHMGSetInnerPCType_HMG(PC pc, PCType type)
347fd2dd295SFande Kong {
34807a4832bSFande Kong   PC_MG           *mg  = (PC_MG*)pc->data;
34907a4832bSFande Kong   PC_HMG          *hmg = (PC_HMG*) mg->innerctx;
350fd2dd295SFande Kong 
351fd2dd295SFande Kong   PetscFunctionBegin;
3529566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(type,&(hmg->innerpctype)));
3538a2c336bSFande Kong   PetscFunctionReturn(0);
3548a2c336bSFande Kong }
3558a2c336bSFande Kong 
356b155ba7fSFande Kong /*@C
3578a2c336bSFande Kong    PCHMGSetInnerPCType - Set an inner PC type
3588a2c336bSFande Kong 
3598a2c336bSFande Kong    Logically Collective on PC
3608a2c336bSFande Kong 
3618a2c336bSFande Kong    Input Parameters:
3628a2c336bSFande Kong +  pc - the HMG context
3638a2c336bSFande Kong -  type - <hypre, gamg> coarsening algorithm
3648a2c336bSFande Kong 
3658a2c336bSFande Kong    Options Database Keys:
36649c604d5SFande Kong .  -hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix
3678a2c336bSFande Kong 
3688a2c336bSFande Kong    Level: beginner
3698a2c336bSFande Kong 
3708a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, coarsening
3718a2c336bSFande Kong 
3728a2c336bSFande Kong .seealso: PCHMG, PCType
373b155ba7fSFande Kong @*/
3748a2c336bSFande Kong PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type)
3758a2c336bSFande Kong {
3768a2c336bSFande Kong   PetscFunctionBegin;
3778a2c336bSFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
378*cac4c232SBarry Smith   PetscUseMethod(pc,"PCHMGSetInnerPCType_C",(PC,PCType),(pc,type));
3798a2c336bSFande Kong   PetscFunctionReturn(0);
3808a2c336bSFande Kong }
3818a2c336bSFande Kong 
38249c604d5SFande Kong static PetscErrorCode PCHMGSetCoarseningComponent_HMG(PC pc, PetscInt component)
38349c604d5SFande Kong {
38449c604d5SFande Kong   PC_MG           *mg  = (PC_MG*)pc->data;
38549c604d5SFande Kong   PC_HMG          *hmg = (PC_HMG*) mg->innerctx;
38649c604d5SFande Kong 
38749c604d5SFande Kong   PetscFunctionBegin;
38849c604d5SFande Kong   hmg->component = component;
38949c604d5SFande Kong   PetscFunctionReturn(0);
39049c604d5SFande Kong }
39149c604d5SFande Kong 
392b155ba7fSFande Kong /*@
39349c604d5SFande Kong    PCHMGSetCoarseningComponent - Set which component is used for the subspace-based coarsening algorithm
39449c604d5SFande Kong 
39549c604d5SFande Kong    Logically Collective on PC
39649c604d5SFande Kong 
39749c604d5SFande Kong    Input Parameters:
39849c604d5SFande Kong +  pc - the HMG context
39949c604d5SFande Kong -  component - which component PC will coarsen
40049c604d5SFande Kong 
40149c604d5SFande Kong    Options Database Keys:
40249c604d5SFande Kong .  -pc_hmg_coarsening_component - Which component is chosen for the subspace-based coarsening algorithm
40349c604d5SFande Kong 
40449c604d5SFande Kong    Level: beginner
40549c604d5SFande Kong 
40649c604d5SFande Kong .keywords: HMG, multigrid, interpolation, coarsening, component
40749c604d5SFande Kong 
40849c604d5SFande Kong .seealso: PCHMG, PCType
409b155ba7fSFande Kong @*/
41049c604d5SFande Kong PetscErrorCode PCHMGSetCoarseningComponent(PC pc, PetscInt component)
41149c604d5SFande Kong {
41249c604d5SFande Kong   PetscFunctionBegin;
41349c604d5SFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
414*cac4c232SBarry Smith   PetscUseMethod(pc,"PCHMGSetCoarseningComponent_C",(PC,PetscInt),(pc,component));
41549c604d5SFande Kong   PetscFunctionReturn(0);
41649c604d5SFande Kong }
41749c604d5SFande Kong 
41849c604d5SFande Kong static PetscErrorCode PCHMGUseMatMAIJ_HMG(PC pc, PetscBool usematmaij)
41949c604d5SFande Kong {
42049c604d5SFande Kong   PC_MG           *mg  = (PC_MG*)pc->data;
42149c604d5SFande Kong   PC_HMG          *hmg = (PC_HMG*) mg->innerctx;
42249c604d5SFande Kong 
42349c604d5SFande Kong   PetscFunctionBegin;
42449c604d5SFande Kong   hmg->usematmaij = usematmaij;
42549c604d5SFande Kong   PetscFunctionReturn(0);
42649c604d5SFande Kong }
42749c604d5SFande Kong 
428b155ba7fSFande Kong /*@
42949c604d5SFande Kong    PCHMGUseMatMAIJ - Set a flag that indicates if or not to use MatMAIJ for interpolations for saving memory
43049c604d5SFande Kong 
43149c604d5SFande Kong    Logically Collective on PC
43249c604d5SFande Kong 
43349c604d5SFande Kong    Input Parameters:
43449c604d5SFande Kong +  pc - the HMG context
43549c604d5SFande Kong -  usematmaij - if or not to use MatMAIJ for interpolations. By default, it is true for saving memory
43649c604d5SFande Kong 
43749c604d5SFande Kong    Options Database Keys:
43849c604d5SFande Kong .  -pc_hmg_use_matmaij - <true | false >
43949c604d5SFande Kong 
44049c604d5SFande Kong    Level: beginner
44149c604d5SFande Kong 
44249c604d5SFande Kong .keywords: HMG, multigrid, interpolation, coarsening, MatMAIJ
44349c604d5SFande Kong 
44449c604d5SFande Kong .seealso: PCHMG, PCType
445b155ba7fSFande Kong @*/
44649c604d5SFande Kong PetscErrorCode PCHMGUseMatMAIJ(PC pc, PetscBool usematmaij)
44749c604d5SFande Kong {
44849c604d5SFande Kong   PetscFunctionBegin;
44949c604d5SFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
450*cac4c232SBarry Smith   PetscUseMethod(pc,"PCHMGUseMatMAIJ_C",(PC,PetscBool),(pc,usematmaij));
45149c604d5SFande Kong   PetscFunctionReturn(0);
45249c604d5SFande Kong }
45349c604d5SFande Kong 
4548a2c336bSFande Kong /*MC
455fd2dd295SFande Kong    PCHMG - Hybrid of PETSc preconditioners (such as ASM, BJacobi, SOR, etc.) and Hypre BoomerAMG, GAMG or other multilevel methods. BoomerAMG, GAMG
456fd2dd295SFande Kong            or other multilevel methods is used to coarsen matrix and generate a sequence of coarse matrices and interpolations. The matrices and
457fd2dd295SFande Kong            interpolations are employed to construct PCMG, and then any available PETSc preconditioners can be chosen as smoothers and the coarse solver.
458360ee056SFande Kong 
459360ee056SFande Kong    Options Database Keys:
4608a2c336bSFande Kong +  -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
4618a2c336bSFande Kong .  -pc_hmg_use_subspace_coarsening  <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix).
462fd2dd295SFande Kong .  -hmg_inner_pc_type <hypre, gamg, ...> - What method is used to coarsen matrix
4634bb91820SFande Kong -  -pc_hmg_use_matmaij <true | false> - Whether or not to use MatMAIJ for multicomponent problems for saving memory
464360ee056SFande Kong 
465360ee056SFande Kong    Notes:
466360ee056SFande Kong     For multicomponent problems, we can just coarsen one submatrix associated with one particular component. In this way, the preconditioner setup
467360ee056SFande Kong     time is significantly reduced. One typical use case is neutron transport equations. There are many variables on each mesh vertex due to the
468360ee056SFande Kong     of angle and energy. Each variable, in fact, corresponds to the same PDEs but with different material properties.
469360ee056SFande Kong 
470360ee056SFande Kong    Level: beginner
471360ee056SFande Kong 
4728a2c336bSFande Kong    Concepts: Hybrid of ASM and MG, Subspace Coarsening
473360ee056SFande Kong 
474360ee056SFande Kong     References:
475606c0280SSatish Balay .   * - Fande Kong, Yaqi Wang, Derek R Gaston, Cody J Permann, Andrew E Slaughter, Alexander D Lindsay, Richard C Martineau, A highly parallel multilevel
476360ee056SFande Kong     Newton-Krylov-Schwarz method with subspace-based coarsening and partition-based balancing for the multigroup neutron transport equations on
477360ee056SFande Kong     3D unstructured meshes, arXiv preprint arXiv:1903.03659, 2019
478360ee056SFande Kong 
479fd2dd295SFande Kong .seealso:  PCCreate(), PCSetType(), PCType, PC, PCMG, PCHYPRE, PCHMG, PCGetCoarseOperators(), PCGetInterpolations(), PCHMGSetReuseInterpolation(), PCHMGSetUseSubspaceCoarsening(),
480fd2dd295SFande Kong            PCHMGSetInnerPCType()
481360ee056SFande Kong 
482360ee056SFande Kong M*/
483360ee056SFande Kong PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc)
484360ee056SFande Kong {
485360ee056SFande Kong   PC_HMG         *hmg;
486360ee056SFande Kong   PC_MG          *mg;
487360ee056SFande Kong 
488360ee056SFande Kong   PetscFunctionBegin;
489360ee056SFande Kong   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
490360ee056SFande Kong   if (pc->ops->destroy) {
4919566063dSJacob Faibussowitsch     PetscCall((*pc->ops->destroy)(pc));
4920a545947SLisandro Dalcin     pc->data = NULL;
493360ee056SFande Kong   }
4949566063dSJacob Faibussowitsch   PetscCall(PetscFree(((PetscObject)pc)->type_name));
495360ee056SFande Kong 
4969566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc,PCMG));
4979566063dSJacob Faibussowitsch   PetscCall(PetscObjectChangeTypeName((PetscObject)pc, PCHMG));
4989566063dSJacob Faibussowitsch   PetscCall(PetscNew(&hmg));
499360ee056SFande Kong 
500360ee056SFande Kong   mg                      = (PC_MG*) pc->data;
501360ee056SFande Kong   mg->innerctx            = hmg;
502360ee056SFande Kong   hmg->reuseinterp        = PETSC_FALSE;
5038a2c336bSFande Kong   hmg->subcoarsening      = PETSC_FALSE;
5044bb91820SFande Kong   hmg->usematmaij         = PETSC_TRUE;
50549c604d5SFande Kong   hmg->component          = 0;
5068a2c336bSFande Kong   hmg->innerpc            = NULL;
507360ee056SFande Kong 
508360ee056SFande Kong   pc->ops->setfromoptions = PCSetFromOptions_HMG;
509360ee056SFande Kong   pc->ops->view           = PCView_HMG;
510360ee056SFande Kong   pc->ops->destroy        = PCDestroy_HMG;
511360ee056SFande Kong   pc->ops->setup          = PCSetUp_HMG;
512fd2dd295SFande Kong 
5139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",PCHMGSetReuseInterpolation_HMG));
5149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",PCHMGSetUseSubspaceCoarsening_HMG));
5159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",PCHMGSetInnerPCType_HMG));
5169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetCoarseningComponent_C",PCHMGSetCoarseningComponent_HMG));
5179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCHMGUseMatMAIJ_C",PCHMGUseMatMAIJ_HMG));
518360ee056SFande Kong   PetscFunctionReturn(0);
519360ee056SFande Kong }
520