xref: /petsc/src/ksp/pc/impls/hmg/hmg.c (revision 49c604d5ad36f286fab770f4b65fe1d26a2afaf8)
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 */
11*49c604d5SFande Kong   PetscBool        subcoarsening;      /* If or not to use a subspace-based coarsening algorithm */
12*49c604d5SFande Kong   PetscBool        usematmaij;         /* If or not to use MatMAIJ for saving memory */
13*49c604d5SFande 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 
19*49c604d5SFande Kong static PetscErrorCode PCHMGExtractSubMatrix_Private(Mat pmat,Mat *submat,MatReuse reuse,PetscInt component,PetscInt blocksize)
20360ee056SFande Kong {
21360ee056SFande Kong   IS             isrow;
22360ee056SFande Kong   PetscErrorCode ierr;
238a2c336bSFande Kong   PetscInt       rstart,rend;
24360ee056SFande Kong   MPI_Comm       comm;
25360ee056SFande Kong 
26360ee056SFande Kong   PetscFunctionBegin;
27360ee056SFande Kong   ierr = PetscObjectGetComm((PetscObject)pmat,&comm);CHKERRQ(ierr);
28*49c604d5SFande Kong   if (component>=blocksize) SETERRQ3(comm,PETSC_ERR_ARG_INCOMP,"Component %d should be less than block size %d \n",component,blocksize);
29360ee056SFande Kong   ierr = MatGetOwnershipRange(pmat,&rstart,&rend);CHKERRQ(ierr);
30360ee056SFande Kong   if ((rend-rstart)%blocksize != 0) SETERRQ3(comm,PETSC_ERR_ARG_INCOMP,"Block size %d is inconsisent for [%d, %d) \n",blocksize,rstart,rend);
31*49c604d5SFande Kong   ierr = ISCreateStride(comm,(rend-rstart)/blocksize,rstart+component,blocksize,&isrow);
32360ee056SFande Kong   ierr = MatCreateSubMatrix(pmat,isrow,isrow,reuse,submat);CHKERRQ(ierr);
33360ee056SFande Kong   ierr = ISDestroy(&isrow);CHKERRQ(ierr);
34360ee056SFande Kong   PetscFunctionReturn(0);
35360ee056SFande Kong }
36360ee056SFande Kong 
378a2c336bSFande Kong static PetscErrorCode PCHMGExpandInterpolation_Private(Mat subinterp, Mat *interp, PetscInt blocksize)
38360ee056SFande Kong {
39360ee056SFande Kong   PetscInt              subrstart,subrend,subrowsize,subcolsize,subcstart,subcend,rowsize,colsize;
408a2c336bSFande Kong   PetscInt              subrow,row,nz,*d_nnz,*o_nnz,i,j,dnz,onz,max_nz,*indices;
418a2c336bSFande Kong   const PetscInt        *idx;
428a2c336bSFande Kong   const PetscScalar     *values;
43360ee056SFande Kong   PetscErrorCode        ierr;
448a2c336bSFande Kong   MPI_Comm              comm;
45360ee056SFande Kong 
46360ee056SFande Kong   PetscFunctionBegin;
478a2c336bSFande Kong   ierr = PetscObjectGetComm((PetscObject)subinterp,&comm);CHKERRQ(ierr);
488a2c336bSFande Kong   ierr = MatGetOwnershipRange(subinterp,&subrstart,&subrend);CHKERRQ(ierr);
498a2c336bSFande Kong   subrowsize = subrend-subrstart;
50360ee056SFande Kong   rowsize = subrowsize*blocksize;
51360ee056SFande Kong   ierr = PetscCalloc2(rowsize,&d_nnz,rowsize,&o_nnz);CHKERRQ(ierr);
528a2c336bSFande Kong   ierr = MatGetOwnershipRangeColumn(subinterp,&subcstart,&subcend);CHKERRQ(ierr);
538a2c336bSFande Kong   subcolsize = subcend - subcstart;
548a2c336bSFande Kong   colsize    = subcolsize*blocksize;
55360ee056SFande Kong   max_nz = 0;
568a2c336bSFande Kong   for (subrow=subrstart;subrow<subrend;subrow++) {
578a2c336bSFande Kong     ierr = MatGetRow(subinterp,subrow,&nz,&idx,NULL);CHKERRQ(ierr);
58360ee056SFande Kong     if (max_nz<nz) max_nz = nz;
59360ee056SFande Kong     dnz = 0; onz = 0;
60360ee056SFande Kong     for (i=0;i<nz;i++) {
618a2c336bSFande Kong       if(idx[i]>=subcstart && idx[i]<subcend) dnz++;
628a2c336bSFande Kong       else onz++;
63360ee056SFande Kong     }
64360ee056SFande Kong     for (i=0;i<blocksize;i++) {
65360ee056SFande Kong       d_nnz[(subrow-subrstart)*blocksize+i] = dnz;
66360ee056SFande Kong       o_nnz[(subrow-subrstart)*blocksize+i] = onz;
67360ee056SFande Kong     }
688a2c336bSFande Kong     ierr = MatRestoreRow(subinterp,subrow,&nz,&idx,NULL);CHKERRQ(ierr);
69360ee056SFande Kong   }
70360ee056SFande Kong   ierr = MatCreateAIJ(comm,rowsize,colsize,PETSC_DETERMINE,PETSC_DETERMINE,0,d_nnz,0,o_nnz,interp);CHKERRQ(ierr);
71360ee056SFande Kong   ierr = MatSetOption(*interp,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
72360ee056SFande Kong   ierr = MatSetOption(*interp,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
73360ee056SFande Kong   ierr = MatSetOption(*interp,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
74360ee056SFande Kong   ierr = MatSetFromOptions(*interp);CHKERRQ(ierr);
75360ee056SFande Kong 
76360ee056SFande Kong   ierr = MatSetUp(*interp);CHKERRQ(ierr);
77360ee056SFande Kong   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
78071fcb05SBarry Smith   ierr = PetscMalloc1(max_nz,&indices);CHKERRQ(ierr);
798a2c336bSFande Kong   for (subrow=subrstart; subrow<subrend; subrow++) {
808a2c336bSFande Kong     ierr = MatGetRow(subinterp,subrow,&nz,&idx,&values);CHKERRQ(ierr);
81360ee056SFande Kong     for (i=0;i<blocksize;i++) {
82360ee056SFande Kong       row = subrow*blocksize+i;
83360ee056SFande Kong       for (j=0;j<nz;j++) {
84360ee056SFande Kong         indices[j] = idx[j]*blocksize+i;
85360ee056SFande Kong       }
86360ee056SFande Kong       ierr = MatSetValues(*interp,1,&row,nz,indices,values,INSERT_VALUES);CHKERRQ(ierr);
87360ee056SFande Kong     }
888a2c336bSFande Kong     ierr = MatRestoreRow(subinterp,subrow,&nz,&idx,&values);CHKERRQ(ierr);
89360ee056SFande Kong   }
90360ee056SFande Kong   ierr = PetscFree(indices);CHKERRQ(ierr);
91360ee056SFande Kong   ierr = MatAssemblyBegin(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92360ee056SFande Kong   ierr = MatAssemblyEnd(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
93360ee056SFande Kong   PetscFunctionReturn(0);
94360ee056SFande Kong }
95360ee056SFande Kong 
96360ee056SFande Kong PetscErrorCode PCSetUp_HMG(PC pc)
97360ee056SFande Kong {
98360ee056SFande Kong   PetscErrorCode     ierr;
99360ee056SFande Kong   Mat                PA, submat;
100360ee056SFande Kong   PC_MG              *mg   = (PC_MG*)pc->data;
101360ee056SFande Kong   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
102360ee056SFande Kong   MPI_Comm           comm;
103360ee056SFande Kong   PetscInt           level;
104360ee056SFande Kong   PetscInt           num_levels;
1058a2c336bSFande Kong   Mat                *operators,*interpolations;
1068a2c336bSFande Kong   PetscInt           blocksize;
107fd2dd295SFande Kong   const char         *prefix;
10807a4832bSFande Kong   PCMGGalerkinType   galerkin;
109360ee056SFande Kong 
110360ee056SFande Kong   PetscFunctionBegin;
111360ee056SFande Kong   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
112360ee056SFande Kong   if (pc->setupcalled) {
1135f36d8f2SFande Kong    if (hmg->reuseinterp) {
1145f36d8f2SFande Kong      /* If we did not use Galerkin in the last call or we have a different sparsity pattern now,
1155f36d8f2SFande Kong       * we have to build from scratch
11607a4832bSFande Kong       * */
11707a4832bSFande Kong      ierr = PCMGGetGalerkin(pc,&galerkin);CHKERRQ(ierr);
1185f36d8f2SFande Kong      if (galerkin == PC_MG_GALERKIN_NONE || pc->flag != SAME_NONZERO_PATTERN) pc->setupcalled = PETSC_FALSE;
119360ee056SFande Kong      ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT);CHKERRQ(ierr);
120360ee056SFande Kong      ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
121360ee056SFande Kong      PetscFunctionReturn(0);
122360ee056SFande Kong     } else {
123360ee056SFande Kong      ierr = PCReset_MG(pc);CHKERRQ(ierr);
124360ee056SFande Kong      pc->setupcalled = PETSC_FALSE;
125360ee056SFande Kong     }
126360ee056SFande Kong   }
127360ee056SFande Kong 
1288a2c336bSFande Kong   /* Create an inner PC (GAMG or HYPRE) */
1298a2c336bSFande Kong   if (!hmg->innerpc) {
1308a2c336bSFande Kong     ierr = PCCreate(comm,&hmg->innerpc);CHKERRQ(ierr);
131fd2dd295SFande Kong     /* If users do not set an inner pc type, we need to set a default value */
132fd2dd295SFande Kong     if (!hmg->innerpctype) {
133fd2dd295SFande Kong     /* If hypre is available, use hypre, otherwise, use gamg */
134fd2dd295SFande Kong #if PETSC_HAVE_HYPRE
135fd2dd295SFande Kong       ierr = PetscStrallocpy(PCHYPRE,&(hmg->innerpctype));CHKERRQ(ierr);
136fd2dd295SFande Kong #else
137fd2dd295SFande Kong       ierr = PetscStrallocpy(PCGAMG,&(hmg->innerpctype));CHKERRQ(ierr);
138fd2dd295SFande Kong #endif
1398a2c336bSFande Kong     }
140fd2dd295SFande Kong     ierr = PCSetType(hmg->innerpc,hmg->innerpctype);CHKERRQ(ierr);
141360ee056SFande Kong   }
142360ee056SFande Kong   ierr = PCGetOperators(pc,NULL,&PA);CHKERRQ(ierr);
1438a2c336bSFande Kong   /* Users need to correctly set a block size of matrix in order to use subspace coarsening */
1448a2c336bSFande Kong   ierr = MatGetBlockSize(PA,&blocksize);CHKERRQ(ierr);
1458a2c336bSFande Kong   if (blocksize<=1) hmg->subcoarsening = PETSC_FALSE;
1468a2c336bSFande Kong   /* Extract a submatrix for constructing subinterpolations */
1478a2c336bSFande Kong   if (hmg->subcoarsening) {
148*49c604d5SFande Kong     ierr = PCHMGExtractSubMatrix_Private(PA,&submat,MAT_INITIAL_MATRIX,hmg->component,blocksize);CHKERRQ(ierr);
149360ee056SFande Kong     PA = submat;
150360ee056SFande Kong   }
1518a2c336bSFande Kong   ierr = PCSetOperators(hmg->innerpc,PA,PA);CHKERRQ(ierr);
1528a2c336bSFande Kong   if (hmg->subcoarsening) {
153360ee056SFande Kong    ierr = MatDestroy(&PA);CHKERRQ(ierr);
154360ee056SFande Kong   }
1558a2c336bSFande Kong   /* Setup inner PC correctly. During this step, matrix will be coarsened */
1568a2c336bSFande Kong   ierr = PCSetUseAmat(hmg->innerpc,PETSC_FALSE);CHKERRQ(ierr);
157fd2dd295SFande Kong   ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,&prefix);CHKERRQ(ierr);
158fd2dd295SFande Kong   ierr = PetscObjectSetOptionsPrefix((PetscObject)hmg->innerpc,prefix);CHKERRQ(ierr);
159fd2dd295SFande Kong   ierr = PetscObjectAppendOptionsPrefix((PetscObject)hmg->innerpc,"hmg_inner_");CHKERRQ(ierr);
160fd2dd295SFande Kong   ierr = PCSetFromOptions(hmg->innerpc);CHKERRQ(ierr);
1618a2c336bSFande Kong   ierr = PCSetUp(hmg->innerpc);
162360ee056SFande Kong 
1638a2c336bSFande Kong   /* Obtain interpolations IN PLACE. For BoomerAMG, (I,J,data) is reused to avoid memory overhead */
164fd2dd295SFande Kong   ierr = PCGetInterpolations(hmg->innerpc,&num_levels,&interpolations);CHKERRQ(ierr);
1658a2c336bSFande Kong   /* We can reuse the coarse operators when we do the full space coarsening */
1668a2c336bSFande Kong   if (!hmg->subcoarsening) {
167fd2dd295SFande Kong     ierr = PCGetCoarseOperators(hmg->innerpc,&num_levels,&operators);CHKERRQ(ierr);
168360ee056SFande Kong   }
169360ee056SFande Kong 
1708a2c336bSFande Kong   ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr);
1718a2c336bSFande Kong   hmg->innerpc = 0;
1728a2c336bSFande Kong   ierr = PCMGSetLevels_MG(pc,num_levels,NULL);CHKERRQ(ierr);
1738a2c336bSFande Kong   /* Set coarse matrices and interpolations to PCMG */
174360ee056SFande Kong   for (level=num_levels-1; level>0; level--) {
175360ee056SFande Kong     Mat P=0, pmat=0;
176360ee056SFande Kong     Vec b, x,r;
1778a2c336bSFande Kong     if (hmg->subcoarsening) {
1784bb91820SFande Kong       if (hmg->usematmaij) {
1794bb91820SFande Kong         ierr = MatCreateMAIJ(interpolations[level-1],blocksize,&P);CHKERRQ(ierr);
1804bb91820SFande Kong         ierr = MatDestroy(&interpolations[level-1]);CHKERRQ(ierr);
1814bb91820SFande Kong       } else {
1828a2c336bSFande Kong         /* Grow interpolation. In the future, we should use MAIJ */
1838a2c336bSFande Kong         ierr = PCHMGExpandInterpolation_Private(interpolations[level-1],&P,blocksize);CHKERRQ(ierr);
1848a2c336bSFande Kong         ierr = MatDestroy(&interpolations[level-1]);CHKERRQ(ierr);
1854bb91820SFande Kong       }
186360ee056SFande Kong     } else {
1878a2c336bSFande Kong       P = interpolations[level-1];
188360ee056SFande Kong     }
189360ee056SFande Kong     ierr = MatCreateVecs(P,&b,&r);CHKERRQ(ierr);
190360ee056SFande Kong     ierr = PCMGSetInterpolation(pc,level,P);CHKERRQ(ierr);
191360ee056SFande Kong     ierr = PCMGSetRestriction(pc,level,P);CHKERRQ(ierr);
192360ee056SFande Kong     ierr = MatDestroy(&P);CHKERRQ(ierr);
1938a2c336bSFande Kong     /* We reuse the matrices when we do not do subspace coarsening */
1948a2c336bSFande Kong     if ((level-1)>=0 && !hmg->subcoarsening) {
1958a2c336bSFande Kong       pmat = operators[level-1];
1968a2c336bSFande Kong       ierr = PCMGSetOperators(pc,level-1,pmat,pmat);CHKERRQ(ierr);
197360ee056SFande Kong       ierr = MatDestroy(&pmat);CHKERRQ(ierr);
198360ee056SFande Kong     }
199360ee056SFande Kong     ierr = PCMGSetRhs(pc,level-1,b);CHKERRQ(ierr);
200360ee056SFande Kong 
201360ee056SFande Kong     ierr = PCMGSetR(pc,level,r);CHKERRQ(ierr);
202360ee056SFande Kong     ierr = VecDestroy(&r);CHKERRQ(ierr);
203360ee056SFande Kong 
204fd2dd295SFande Kong     ierr = VecDuplicate(b,&x);CHKERRQ(ierr);
205360ee056SFande Kong     ierr = PCMGSetX(pc,level-1,x);CHKERRQ(ierr);
206360ee056SFande Kong     ierr = VecDestroy(&x);CHKERRQ(ierr);
207360ee056SFande Kong     ierr = VecDestroy(&b);CHKERRQ(ierr);
208360ee056SFande Kong   }
2098a2c336bSFande Kong   ierr = PetscFree(interpolations);CHKERRQ(ierr);
2108a2c336bSFande Kong   if (!hmg->subcoarsening) {
2118a2c336bSFande Kong     ierr = PetscFree(operators);CHKERRQ(ierr);
2128a2c336bSFande Kong   }
2138a2c336bSFande Kong   /* Turn Galerkin off when we already have coarse operators */
2148a2c336bSFande Kong   ierr = PCMGSetGalerkin(pc,hmg->subcoarsening ? PC_MG_GALERKIN_PMAT:PC_MG_GALERKIN_NONE);CHKERRQ(ierr);
215360ee056SFande Kong   ierr = PCSetDM(pc,NULL);CHKERRQ(ierr);
216360ee056SFande Kong   ierr = PCSetUseAmat(pc,PETSC_FALSE);CHKERRQ(ierr);
217360ee056SFande Kong   ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
218360ee056SFande Kong   ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */
219360ee056SFande Kong   ierr = PetscOptionsEnd();CHKERRQ(ierr);
220360ee056SFande Kong   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
221360ee056SFande Kong   PetscFunctionReturn(0);
222360ee056SFande Kong }
223360ee056SFande Kong 
224360ee056SFande Kong PetscErrorCode PCDestroy_HMG(PC pc)
225360ee056SFande Kong {
226360ee056SFande Kong   PetscErrorCode ierr;
227360ee056SFande Kong   PC_MG          *mg  = (PC_MG*)pc->data;
2288a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
229360ee056SFande Kong 
230360ee056SFande Kong   PetscFunctionBegin;
2318a2c336bSFande Kong   ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr);
2328a2c336bSFande Kong   ierr = PetscFree(hmg->innerpctype);CHKERRQ(ierr);
2338a2c336bSFande Kong   ierr = PetscFree(hmg);CHKERRQ(ierr);
234360ee056SFande Kong   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
235fd2dd295SFande Kong 
236fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",NULL);CHKERRQ(ierr);
237fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",NULL);CHKERRQ(ierr);
238fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",NULL);CHKERRQ(ierr);
239360ee056SFande Kong   PetscFunctionReturn(0);
240360ee056SFande Kong }
241360ee056SFande Kong 
242360ee056SFande Kong PetscErrorCode PCView_HMG(PC pc,PetscViewer viewer)
243360ee056SFande Kong {
244360ee056SFande Kong   PC_MG          *mg = (PC_MG*)pc->data;
2458a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
246360ee056SFande Kong   PetscErrorCode ierr;
247360ee056SFande Kong   PetscBool      iascii;
248360ee056SFande Kong 
249360ee056SFande Kong   PetscFunctionBegin;
250360ee056SFande Kong   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
251360ee056SFande Kong   if (iascii) {
2528a2c336bSFande Kong     ierr = PetscViewerASCIIPrintf(viewer," Reuse interpolation: %s\n",hmg->reuseinterp? "true":"false");CHKERRQ(ierr);
2538a2c336bSFande Kong     ierr = PetscViewerASCIIPrintf(viewer," Use subspace coarsening: %s\n",hmg->subcoarsening? "true":"false");CHKERRQ(ierr);
2548a2c336bSFande Kong     ierr = PetscViewerASCIIPrintf(viewer," Inner PC type: %s \n",hmg->innerpctype);CHKERRQ(ierr);
255360ee056SFande Kong   }
256360ee056SFande Kong   ierr = PCView_MG(pc,viewer);CHKERRQ(ierr);
257360ee056SFande Kong   PetscFunctionReturn(0);
258360ee056SFande Kong }
259360ee056SFande Kong 
260360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems *PetscOptionsObject,PC pc)
261360ee056SFande Kong {
262360ee056SFande Kong   PetscErrorCode ierr;
263360ee056SFande Kong   PC_MG          *mg = (PC_MG*)pc->data;
2648a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
265360ee056SFande Kong 
266360ee056SFande Kong   PetscFunctionBegin;
267360ee056SFande Kong   ierr = PetscOptionsHead(PetscOptionsObject,"HMG");CHKERRQ(ierr);
2688a2c336bSFande Kong   ierr = PetscOptionsBool("-pc_hmg_reuse_interpolation","Reuse the interpolation operators when possible (cheaper, weaker when matrix entries change a lot)","None",hmg->reuseinterp,&hmg->reuseinterp,NULL);CHKERRQ(ierr);
2698a2c336bSFande Kong   ierr = PetscOptionsBool("-pc_hmg_use_subspace_coarsening","Use the subspace coarsening to compute the interpolations","None",hmg->subcoarsening,&hmg->subcoarsening,NULL);CHKERRQ(ierr);
2704bb91820SFande Kong   ierr = PetscOptionsBool("-pc_hmg_use_matmaij","Use MatMAIJ store interpolation for saving memory","None",hmg->usematmaij,&hmg->usematmaij,NULL);CHKERRQ(ierr);
271*49c604d5SFande Kong   ierr = PetscOptionsBool("-pc_hmg_coarsening_component","Which component is chosen for the subspace-based coarsening algorithm","None",hmg->component,&hmg->component,NULL);CHKERRQ(ierr);
272360ee056SFande Kong   ierr = PetscOptionsTail();CHKERRQ(ierr);
273360ee056SFande Kong   PetscFunctionReturn(0);
274360ee056SFande Kong }
275360ee056SFande Kong 
276fd2dd295SFande Kong static PetscErrorCode PCHMGSetReuseInterpolation_HMG(PC pc, PetscBool reuse)
277fd2dd295SFande Kong {
27807a4832bSFande Kong   PC_MG          *mg  = (PC_MG*)pc->data;
27907a4832bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
280fd2dd295SFande Kong 
281fd2dd295SFande Kong   PetscFunctionBegin;
282fd2dd295SFande Kong   hmg->reuseinterp = reuse;
283fd2dd295SFande Kong   PetscFunctionReturn(0);
284fd2dd295SFande Kong }
285fd2dd295SFande Kong 
286360ee056SFande Kong /*MC
2878a2c336bSFande Kong    PCHMGSetReuseInterpolation - Reuse interpolation matrices in HMG
2888a2c336bSFande Kong 
2898a2c336bSFande Kong    Logically Collective on PC
2908a2c336bSFande Kong 
2918a2c336bSFande Kong    Input Parameters:
2928a2c336bSFande Kong +  pc - the HMG context
2938a2c336bSFande Kong -  reuse - True indicates that HMG will reuse the interpolations
2948a2c336bSFande Kong 
2958a2c336bSFande Kong    Options Database Keys:
296*49c604d5SFande Kong .  -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
2978a2c336bSFande Kong 
2988a2c336bSFande Kong    Level: beginner
2998a2c336bSFande Kong 
3008a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, reuse, set
3018a2c336bSFande Kong 
3028a2c336bSFande Kong .seealso: PCHMG
3038a2c336bSFande Kong M*/
3048a2c336bSFande Kong PetscErrorCode PCHMGSetReuseInterpolation(PC pc, PetscBool reuse)
3058a2c336bSFande Kong {
306fd2dd295SFande Kong   PetscErrorCode ierr;
307fd2dd295SFande Kong 
308fd2dd295SFande Kong   PetscFunctionBegin;
309fd2dd295SFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
310fd2dd295SFande Kong   ierr = PetscUseMethod(pc,"PCHMGSetReuseInterpolation_C",(PC,PetscBool),(pc,reuse));CHKERRQ(ierr);
311fd2dd295SFande Kong   PetscFunctionReturn(0);
312fd2dd295SFande Kong }
313fd2dd295SFande Kong 
314fd2dd295SFande Kong static PetscErrorCode PCHMGSetUseSubspaceCoarsening_HMG(PC pc, PetscBool subspace)
315fd2dd295SFande Kong {
31607a4832bSFande Kong   PC_MG          *mg  = (PC_MG*)pc->data;
31707a4832bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
3188a2c336bSFande Kong 
3198a2c336bSFande Kong   PetscFunctionBegin;
320fd2dd295SFande Kong   hmg->subcoarsening = subspace;
3218a2c336bSFande Kong   PetscFunctionReturn(0);
3228a2c336bSFande Kong }
3238a2c336bSFande Kong 
3248a2c336bSFande Kong /*MC
3258a2c336bSFande Kong    PCHMGSetUseSubspaceCoarsening - Use subspace coarsening in HMG
3268a2c336bSFande Kong 
3278a2c336bSFande Kong    Logically Collective on PC
3288a2c336bSFande Kong 
3298a2c336bSFande Kong    Input Parameters:
3308a2c336bSFande Kong +  pc - the HMG context
3318a2c336bSFande Kong -  reuse - True indicates that HMG will use the subspace coarsening
3328a2c336bSFande Kong 
3338a2c336bSFande Kong    Options Database Keys:
334*49c604d5SFande Kong .  -pc_hmg_use_subspace_coarsening  <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix).
3358a2c336bSFande Kong 
3368a2c336bSFande Kong    Level: beginner
3378a2c336bSFande Kong 
3388a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, subspace, coarsening
3398a2c336bSFande Kong 
3408a2c336bSFande Kong .seealso: PCHMG
3418a2c336bSFande Kong M*/
3428a2c336bSFande Kong PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC pc, PetscBool subspace)
3438a2c336bSFande Kong {
344fd2dd295SFande Kong   PetscErrorCode ierr;
3458a2c336bSFande Kong 
3468a2c336bSFande Kong   PetscFunctionBegin;
3478a2c336bSFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
348fd2dd295SFande Kong   ierr = PetscUseMethod(pc,"PCHMGSetUseSubspaceCoarsening_C",(PC,PetscBool),(pc,subspace));CHKERRQ(ierr);
349fd2dd295SFande Kong   PetscFunctionReturn(0);
350fd2dd295SFande Kong }
351fd2dd295SFande Kong 
352fd2dd295SFande Kong static PetscErrorCode PCHMGSetInnerPCType_HMG(PC pc, PCType type)
353fd2dd295SFande Kong {
35407a4832bSFande Kong   PC_MG           *mg  = (PC_MG*)pc->data;
35507a4832bSFande Kong   PC_HMG          *hmg = (PC_HMG*) mg->innerctx;
356fd2dd295SFande Kong   PetscErrorCode  ierr;
357fd2dd295SFande Kong 
358fd2dd295SFande Kong   PetscFunctionBegin;
359fd2dd295SFande Kong   ierr = PetscStrallocpy(type,&(hmg->innerpctype));CHKERRQ(ierr);
3608a2c336bSFande Kong   PetscFunctionReturn(0);
3618a2c336bSFande Kong }
3628a2c336bSFande Kong 
3638a2c336bSFande Kong /*MC
3648a2c336bSFande Kong    PCHMGSetInnerPCType - Set an inner PC type
3658a2c336bSFande Kong 
3668a2c336bSFande Kong    Logically Collective on PC
3678a2c336bSFande Kong 
3688a2c336bSFande Kong    Input Parameters:
3698a2c336bSFande Kong +  pc - the HMG context
3708a2c336bSFande Kong -  type - <hypre, gamg> coarsening algorithm
3718a2c336bSFande Kong 
3728a2c336bSFande Kong    Options Database Keys:
373*49c604d5SFande Kong .  -hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix
3748a2c336bSFande Kong 
3758a2c336bSFande Kong    Level: beginner
3768a2c336bSFande Kong 
3778a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, coarsening
3788a2c336bSFande Kong 
3798a2c336bSFande Kong .seealso: PCHMG, PCType
3808a2c336bSFande Kong M*/
3818a2c336bSFande Kong PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type)
3828a2c336bSFande Kong {
3838a2c336bSFande Kong   PetscErrorCode  ierr;
3848a2c336bSFande Kong 
3858a2c336bSFande Kong   PetscFunctionBegin;
3868a2c336bSFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
387fd2dd295SFande Kong   ierr = PetscUseMethod(pc,"PCHMGSetInnerPCType_C",(PC,PCType),(pc,type));CHKERRQ(ierr);
3888a2c336bSFande Kong   PetscFunctionReturn(0);
3898a2c336bSFande Kong }
3908a2c336bSFande Kong 
391*49c604d5SFande Kong static PetscErrorCode PCHMGSetCoarseningComponent_HMG(PC pc, PetscInt component)
392*49c604d5SFande Kong {
393*49c604d5SFande Kong   PC_MG           *mg  = (PC_MG*)pc->data;
394*49c604d5SFande Kong   PC_HMG          *hmg = (PC_HMG*) mg->innerctx;
395*49c604d5SFande Kong 
396*49c604d5SFande Kong   PetscFunctionBegin;
397*49c604d5SFande Kong   hmg->component = component;
398*49c604d5SFande Kong   PetscFunctionReturn(0);
399*49c604d5SFande Kong }
400*49c604d5SFande Kong 
401*49c604d5SFande Kong /*MC
402*49c604d5SFande Kong    PCHMGSetCoarseningComponent - Set which component is used for the subspace-based coarsening algorithm
403*49c604d5SFande Kong 
404*49c604d5SFande Kong    Logically Collective on PC
405*49c604d5SFande Kong 
406*49c604d5SFande Kong    Input Parameters:
407*49c604d5SFande Kong +  pc - the HMG context
408*49c604d5SFande Kong -  component - which component PC will coarsen
409*49c604d5SFande Kong 
410*49c604d5SFande Kong    Options Database Keys:
411*49c604d5SFande Kong .  -pc_hmg_coarsening_component - Which component is chosen for the subspace-based coarsening algorithm
412*49c604d5SFande Kong 
413*49c604d5SFande Kong    Level: beginner
414*49c604d5SFande Kong 
415*49c604d5SFande Kong .keywords: HMG, multigrid, interpolation, coarsening, component
416*49c604d5SFande Kong 
417*49c604d5SFande Kong .seealso: PCHMG, PCType
418*49c604d5SFande Kong M*/
419*49c604d5SFande Kong PetscErrorCode PCHMGSetCoarseningComponent(PC pc, PetscInt component)
420*49c604d5SFande Kong {
421*49c604d5SFande Kong   PetscErrorCode  ierr;
422*49c604d5SFande Kong 
423*49c604d5SFande Kong   PetscFunctionBegin;
424*49c604d5SFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
425*49c604d5SFande Kong   ierr = PetscUseMethod(pc,"PCHMGSetCoarseningComponent_C",(PC,PetscInt),(pc,component));CHKERRQ(ierr);
426*49c604d5SFande Kong   PetscFunctionReturn(0);
427*49c604d5SFande Kong }
428*49c604d5SFande Kong 
429*49c604d5SFande Kong static PetscErrorCode PCHMGUseMatMAIJ_HMG(PC pc, PetscBool usematmaij)
430*49c604d5SFande Kong {
431*49c604d5SFande Kong   PC_MG           *mg  = (PC_MG*)pc->data;
432*49c604d5SFande Kong   PC_HMG          *hmg = (PC_HMG*) mg->innerctx;
433*49c604d5SFande Kong 
434*49c604d5SFande Kong   PetscFunctionBegin;
435*49c604d5SFande Kong   hmg->usematmaij = usematmaij;
436*49c604d5SFande Kong   PetscFunctionReturn(0);
437*49c604d5SFande Kong }
438*49c604d5SFande Kong 
439*49c604d5SFande Kong /*MC
440*49c604d5SFande Kong    PCHMGUseMatMAIJ - Set a flag that indicates if or not to use MatMAIJ for interpolations for saving memory
441*49c604d5SFande Kong 
442*49c604d5SFande Kong    Logically Collective on PC
443*49c604d5SFande Kong 
444*49c604d5SFande Kong    Input Parameters:
445*49c604d5SFande Kong +  pc - the HMG context
446*49c604d5SFande Kong -  usematmaij - if or not to use MatMAIJ for interpolations. By default, it is true for saving memory
447*49c604d5SFande Kong 
448*49c604d5SFande Kong    Options Database Keys:
449*49c604d5SFande Kong .  -pc_hmg_use_matmaij - <true | false >
450*49c604d5SFande Kong 
451*49c604d5SFande Kong    Level: beginner
452*49c604d5SFande Kong 
453*49c604d5SFande Kong .keywords: HMG, multigrid, interpolation, coarsening, MatMAIJ
454*49c604d5SFande Kong 
455*49c604d5SFande Kong .seealso: PCHMG, PCType
456*49c604d5SFande Kong M*/
457*49c604d5SFande Kong PetscErrorCode PCHMGUseMatMAIJ(PC pc, PetscBool usematmaij)
458*49c604d5SFande Kong {
459*49c604d5SFande Kong   PetscErrorCode  ierr;
460*49c604d5SFande Kong 
461*49c604d5SFande Kong   PetscFunctionBegin;
462*49c604d5SFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
463*49c604d5SFande Kong   ierr = PetscUseMethod(pc,"PCHMGUseMatMAIJ_C",(PC,PetscBool),(pc,usematmaij));CHKERRQ(ierr);
464*49c604d5SFande Kong   PetscFunctionReturn(0);
465*49c604d5SFande Kong }
466*49c604d5SFande Kong 
4678a2c336bSFande Kong /*MC
468fd2dd295SFande Kong    PCHMG - Hybrid of PETSc preconditioners (such as ASM, BJacobi, SOR, etc.) and Hypre BoomerAMG, GAMG or other multilevel methods. BoomerAMG, GAMG
469fd2dd295SFande Kong            or other multilevel methods is used to coarsen matrix and generate a sequence of coarse matrices and interpolations. The matrices and
470fd2dd295SFande Kong            interpolations are employed to construct PCMG, and then any available PETSc preconditioners can be chosen as smoothers and the coarse solver.
471360ee056SFande Kong 
472360ee056SFande Kong    Options Database Keys:
4738a2c336bSFande Kong +  -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
4748a2c336bSFande Kong .  -pc_hmg_use_subspace_coarsening  <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix).
475fd2dd295SFande Kong .  -hmg_inner_pc_type <hypre, gamg, ...> - What method is used to coarsen matrix
4764bb91820SFande Kong -  -pc_hmg_use_matmaij <true | false> - Whether or not to use MatMAIJ for multicomponent problems for saving memory
477360ee056SFande Kong 
478360ee056SFande Kong 
479360ee056SFande Kong    Notes:
480360ee056SFande Kong     For multicomponent problems, we can just coarsen one submatrix associated with one particular component. In this way, the preconditioner setup
481360ee056SFande Kong     time is significantly reduced. One typical use case is neutron transport equations. There are many variables on each mesh vertex due to the
482360ee056SFande Kong     of angle and energy. Each variable, in fact, corresponds to the same PDEs but with different material properties.
483360ee056SFande Kong 
484360ee056SFande Kong    Level: beginner
485360ee056SFande Kong 
4868a2c336bSFande Kong    Concepts: Hybrid of ASM and MG, Subspace Coarsening
487360ee056SFande Kong 
488360ee056SFande Kong     References:
489*49c604d5SFande Kong .   1. - Fande Kong, Yaqi Wang, Derek R Gaston, Cody J Permann, Andrew E Slaughter, Alexander D Lindsay, Richard C Martineau, A highly parallel multilevel
490360ee056SFande Kong     Newton-Krylov-Schwarz method with subspace-based coarsening and partition-based balancing for the multigroup neutron transport equations on
491360ee056SFande Kong     3D unstructured meshes, arXiv preprint arXiv:1903.03659, 2019
492360ee056SFande Kong 
493fd2dd295SFande Kong .seealso:  PCCreate(), PCSetType(), PCType, PC, PCMG, PCHYPRE, PCHMG, PCGetCoarseOperators(), PCGetInterpolations(), PCHMGSetReuseInterpolation(), PCHMGSetUseSubspaceCoarsening(),
494fd2dd295SFande Kong            PCHMGSetInnerPCType()
495360ee056SFande Kong 
496360ee056SFande Kong M*/
497360ee056SFande Kong PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc)
498360ee056SFande Kong {
499360ee056SFande Kong   PetscErrorCode ierr;
500360ee056SFande Kong   PC_HMG         *hmg;
501360ee056SFande Kong   PC_MG          *mg;
502360ee056SFande Kong 
503360ee056SFande Kong   PetscFunctionBegin;
504360ee056SFande Kong   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
505360ee056SFande Kong   if (pc->ops->destroy) {
506360ee056SFande Kong     ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr);
507360ee056SFande Kong     pc->data = 0;
508360ee056SFande Kong   }
509360ee056SFande Kong   ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr);
510360ee056SFande Kong 
511360ee056SFande Kong   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr);
512645b8336SFande Kong   ierr = PetscObjectChangeTypeName((PetscObject)pc, PCHMG);CHKERRQ(ierr);
513fd2dd295SFande Kong   ierr = PetscNew(&hmg);CHKERRQ(ierr);
514360ee056SFande Kong 
515360ee056SFande Kong   mg                      = (PC_MG*) pc->data;
516360ee056SFande Kong   mg->innerctx            = hmg;
517360ee056SFande Kong   hmg->reuseinterp        = PETSC_FALSE;
5188a2c336bSFande Kong   hmg->subcoarsening      = PETSC_FALSE;
5194bb91820SFande Kong   hmg->usematmaij         = PETSC_TRUE;
520*49c604d5SFande Kong   hmg->component          = 0;
5218a2c336bSFande Kong   hmg->innerpc            = NULL;
522360ee056SFande Kong 
523360ee056SFande Kong   pc->ops->setfromoptions = PCSetFromOptions_HMG;
524360ee056SFande Kong   pc->ops->view           = PCView_HMG;
525360ee056SFande Kong   pc->ops->destroy        = PCDestroy_HMG;
526360ee056SFande Kong   pc->ops->setup          = PCSetUp_HMG;
527fd2dd295SFande Kong 
528fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",PCHMGSetReuseInterpolation_HMG);CHKERRQ(ierr);
529fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",PCHMGSetUseSubspaceCoarsening_HMG);CHKERRQ(ierr);
530fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",PCHMGSetInnerPCType_HMG);CHKERRQ(ierr);
531*49c604d5SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetCoarseningComponent_C",PCHMGSetCoarseningComponent_HMG);CHKERRQ(ierr);
532*49c604d5SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGUseMatMAIJ_C",PCHMGUseMatMAIJ_HMG);CHKERRQ(ierr);
533360ee056SFande Kong   PetscFunctionReturn(0);
534360ee056SFande Kong }
535