xref: /petsc/src/ksp/pc/impls/hmg/hmg.c (revision 07a4832b83bce8b7163a7b18eee78858c8f61ab5)
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 */
118a2c336bSFande Kong   PetscBool        subcoarsening;
12360ee056SFande Kong } PC_HMG;
13360ee056SFande Kong 
14360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems*,PC);
15360ee056SFande Kong PetscErrorCode PCReset_MG(PC);
168a2c336bSFande Kong 
178a2c336bSFande Kong static PetscErrorCode PCHMGExtractSubMatrix_Private(Mat pmat,Mat *submat,MatReuse reuse,PetscInt blocksize)
18360ee056SFande Kong {
19360ee056SFande Kong   IS             isrow;
20360ee056SFande Kong   PetscErrorCode ierr;
218a2c336bSFande Kong   PetscInt       rstart,rend;
22360ee056SFande Kong   MPI_Comm       comm;
23360ee056SFande Kong 
24360ee056SFande Kong   PetscFunctionBegin;
25360ee056SFande Kong   ierr = PetscObjectGetComm((PetscObject)pmat,&comm);CHKERRQ(ierr);
26360ee056SFande Kong   ierr = MatGetOwnershipRange(pmat,&rstart,&rend);CHKERRQ(ierr);
27360ee056SFande Kong   if ((rend-rstart)%blocksize != 0) SETERRQ3(comm,PETSC_ERR_ARG_INCOMP,"Block size %d is inconsisent for [%d, %d) \n",blocksize,rstart,rend);
288a2c336bSFande Kong   ierr = ISCreateStride(comm,(rend-rstart)/blocksize,rstart,blocksize,&isrow);
29360ee056SFande Kong   ierr = MatCreateSubMatrix(pmat,isrow,isrow,reuse,submat);CHKERRQ(ierr);
30360ee056SFande Kong   ierr = ISDestroy(&isrow);CHKERRQ(ierr);
31360ee056SFande Kong   PetscFunctionReturn(0);
32360ee056SFande Kong }
33360ee056SFande Kong 
348a2c336bSFande Kong static PetscErrorCode PCHMGExpandInterpolation_Private(Mat subinterp, Mat *interp, PetscInt blocksize)
35360ee056SFande Kong {
36360ee056SFande Kong   PetscInt              subrstart,subrend,subrowsize,subcolsize,subcstart,subcend,rowsize,colsize;
378a2c336bSFande Kong   PetscInt              subrow,row,nz,*d_nnz,*o_nnz,i,j,dnz,onz,max_nz,*indices;
388a2c336bSFande Kong   const PetscInt        *idx;
398a2c336bSFande Kong   const PetscScalar     *values;
40360ee056SFande Kong   PetscErrorCode        ierr;
418a2c336bSFande Kong   MPI_Comm              comm;
42360ee056SFande Kong 
43360ee056SFande Kong   PetscFunctionBegin;
448a2c336bSFande Kong   ierr = PetscObjectGetComm((PetscObject)subinterp,&comm);CHKERRQ(ierr);
458a2c336bSFande Kong   ierr = MatGetOwnershipRange(subinterp,&subrstart,&subrend);CHKERRQ(ierr);
468a2c336bSFande Kong   subrowsize = subrend-subrstart;
47360ee056SFande Kong   rowsize = subrowsize*blocksize;
48360ee056SFande Kong   ierr = PetscCalloc2(rowsize,&d_nnz,rowsize,&o_nnz);CHKERRQ(ierr);
498a2c336bSFande Kong   ierr = MatGetOwnershipRangeColumn(subinterp,&subcstart,&subcend);CHKERRQ(ierr);
508a2c336bSFande Kong   subcolsize = subcend - subcstart;
518a2c336bSFande Kong   colsize    = subcolsize*blocksize;
52360ee056SFande Kong   max_nz = 0;
538a2c336bSFande Kong   for (subrow=subrstart;subrow<subrend;subrow++) {
548a2c336bSFande Kong     ierr = MatGetRow(subinterp,subrow,&nz,&idx,NULL);CHKERRQ(ierr);
55360ee056SFande Kong     if (max_nz<nz) max_nz = nz;
56360ee056SFande Kong     dnz = 0; onz = 0;
57360ee056SFande Kong     for (i=0;i<nz;i++) {
588a2c336bSFande Kong       if(idx[i]>=subcstart && idx[i]<subcend) dnz++;
598a2c336bSFande Kong       else onz++;
60360ee056SFande Kong     }
61360ee056SFande Kong     for (i=0;i<blocksize;i++) {
62360ee056SFande Kong       d_nnz[(subrow-subrstart)*blocksize+i] = dnz;
63360ee056SFande Kong       o_nnz[(subrow-subrstart)*blocksize+i] = onz;
64360ee056SFande Kong     }
658a2c336bSFande Kong     ierr = MatRestoreRow(subinterp,subrow,&nz,&idx,NULL);CHKERRQ(ierr);
66360ee056SFande Kong   }
67360ee056SFande Kong   ierr = MatCreateAIJ(comm,rowsize,colsize,PETSC_DETERMINE,PETSC_DETERMINE,0,d_nnz,0,o_nnz,interp);CHKERRQ(ierr);
68360ee056SFande Kong   ierr = MatSetOption(*interp,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
69360ee056SFande Kong   ierr = MatSetOption(*interp,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
70360ee056SFande Kong   ierr = MatSetOption(*interp,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
71360ee056SFande Kong   ierr = MatSetFromOptions(*interp);CHKERRQ(ierr);
72360ee056SFande Kong 
73360ee056SFande Kong   ierr = MatSetUp(*interp);CHKERRQ(ierr);
74360ee056SFande Kong   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
75360ee056SFande Kong   ierr = PetscCalloc1(max_nz,&indices);CHKERRQ(ierr);
768a2c336bSFande Kong   for (subrow=subrstart; subrow<subrend; subrow++) {
778a2c336bSFande Kong     ierr = MatGetRow(subinterp,subrow,&nz,&idx,&values);CHKERRQ(ierr);
78360ee056SFande Kong     for (i=0;i<blocksize;i++) {
79360ee056SFande Kong       row = subrow*blocksize+i;
80360ee056SFande Kong       for (j=0;j<nz;j++) {
81360ee056SFande Kong         indices[j] = idx[j]*blocksize+i;
82360ee056SFande Kong       }
83360ee056SFande Kong       ierr = MatSetValues(*interp,1,&row,nz,indices,values,INSERT_VALUES);CHKERRQ(ierr);
84360ee056SFande Kong     }
858a2c336bSFande Kong     ierr = MatRestoreRow(subinterp,subrow,&nz,&idx,&values);CHKERRQ(ierr);
86360ee056SFande Kong   }
87360ee056SFande Kong   ierr = PetscFree(indices);CHKERRQ(ierr);
88360ee056SFande Kong   ierr = MatAssemblyBegin(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89360ee056SFande Kong   ierr = MatAssemblyEnd(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90360ee056SFande Kong   PetscFunctionReturn(0);
91360ee056SFande Kong }
92360ee056SFande Kong 
93360ee056SFande Kong PetscErrorCode PCSetUp_HMG(PC pc)
94360ee056SFande Kong {
95360ee056SFande Kong   PetscErrorCode     ierr;
96360ee056SFande Kong   Mat                PA, submat;
97360ee056SFande Kong   PC_MG              *mg   = (PC_MG*)pc->data;
98360ee056SFande Kong   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
99360ee056SFande Kong   MPI_Comm           comm;
100360ee056SFande Kong   PetscInt           level;
101360ee056SFande Kong   PetscInt           num_levels;
1028a2c336bSFande Kong   Mat                *operators,*interpolations;
1038a2c336bSFande Kong   PetscInt           blocksize;
104fd2dd295SFande Kong   const char         *prefix;
105*07a4832bSFande Kong   PCMGGalerkinType   galerkin;
106360ee056SFande Kong 
107360ee056SFande Kong   PetscFunctionBegin;
108360ee056SFande Kong   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
109360ee056SFande Kong   if (pc->setupcalled) {
110*07a4832bSFande Kong    /* Only reuse interpolations when nonzero pattern does not change
111*07a4832bSFande Kong     * since coarse matrices are already allocated.
112*07a4832bSFande Kong     * */
113360ee056SFande Kong    if (pc->flag == SAME_NONZERO_PATTERN && hmg->reuseinterp) {
114*07a4832bSFande Kong      /* If we did not use Galerkin in the last call, we have to build from scratch */
115*07a4832bSFande Kong      ierr = PCMGGetGalerkin(pc,&galerkin);CHKERRQ(ierr);
116*07a4832bSFande Kong      if (galerkin == PC_MG_GALERKIN_NONE) pc->setupcalled = PETSC_FALSE;
117360ee056SFande Kong      ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT);CHKERRQ(ierr);
118360ee056SFande Kong      ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
119360ee056SFande Kong      PetscFunctionReturn(0);
120360ee056SFande Kong     } else {
121360ee056SFande Kong      ierr = PCReset_MG(pc);CHKERRQ(ierr);
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) {
1288a2c336bSFande Kong     ierr = PCCreate(comm,&hmg->innerpc);CHKERRQ(ierr);
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
133fd2dd295SFande Kong       ierr = PetscStrallocpy(PCHYPRE,&(hmg->innerpctype));CHKERRQ(ierr);
134fd2dd295SFande Kong #else
135fd2dd295SFande Kong       ierr = PetscStrallocpy(PCGAMG,&(hmg->innerpctype));CHKERRQ(ierr);
136fd2dd295SFande Kong #endif
1378a2c336bSFande Kong     }
138fd2dd295SFande Kong     ierr = PCSetType(hmg->innerpc,hmg->innerpctype);CHKERRQ(ierr);
139360ee056SFande Kong   }
140360ee056SFande Kong   ierr = PCGetOperators(pc,NULL,&PA);CHKERRQ(ierr);
1418a2c336bSFande Kong   /* Users need to correctly set a block size of matrix in order to use subspace coarsening */
1428a2c336bSFande Kong   ierr = MatGetBlockSize(PA,&blocksize);CHKERRQ(ierr);
1438a2c336bSFande Kong   if (blocksize<=1) hmg->subcoarsening = PETSC_FALSE;
1448a2c336bSFande Kong   /* Extract a submatrix for constructing subinterpolations */
1458a2c336bSFande Kong   if (hmg->subcoarsening) {
1468a2c336bSFande Kong     ierr = PCHMGExtractSubMatrix_Private(PA,&submat,MAT_INITIAL_MATRIX,blocksize);CHKERRQ(ierr);
147360ee056SFande Kong     PA = submat;
148360ee056SFande Kong   }
1498a2c336bSFande Kong   ierr = PCSetOperators(hmg->innerpc,PA,PA);CHKERRQ(ierr);
1508a2c336bSFande Kong   if (hmg->subcoarsening) {
151360ee056SFande Kong    ierr = MatDestroy(&PA);CHKERRQ(ierr);
152360ee056SFande Kong   }
1538a2c336bSFande Kong   /* Setup inner PC correctly. During this step, matrix will be coarsened */
1548a2c336bSFande Kong   ierr = PCSetUseAmat(hmg->innerpc,PETSC_FALSE);CHKERRQ(ierr);
155fd2dd295SFande Kong   ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,&prefix);CHKERRQ(ierr);
156fd2dd295SFande Kong   ierr = PetscObjectSetOptionsPrefix((PetscObject)hmg->innerpc,prefix);CHKERRQ(ierr);
157fd2dd295SFande Kong   ierr = PetscObjectAppendOptionsPrefix((PetscObject)hmg->innerpc,"hmg_inner_");CHKERRQ(ierr);
158fd2dd295SFande Kong   ierr = PCSetFromOptions(hmg->innerpc);CHKERRQ(ierr);
1598a2c336bSFande Kong   ierr = PCSetUp(hmg->innerpc);
160360ee056SFande Kong 
1618a2c336bSFande Kong   /* Obtain interpolations IN PLACE. For BoomerAMG, (I,J,data) is reused to avoid memory overhead */
162fd2dd295SFande Kong   ierr = PCGetInterpolations(hmg->innerpc,&num_levels,&interpolations);CHKERRQ(ierr);
1638a2c336bSFande Kong   /* We can reuse the coarse operators when we do the full space coarsening */
1648a2c336bSFande Kong   if (!hmg->subcoarsening) {
165fd2dd295SFande Kong     ierr = PCGetCoarseOperators(hmg->innerpc,&num_levels,&operators);CHKERRQ(ierr);
166360ee056SFande Kong   }
167360ee056SFande Kong 
1688a2c336bSFande Kong   ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr);
1698a2c336bSFande Kong   hmg->innerpc = 0;
1708a2c336bSFande Kong   ierr = PCMGSetLevels_MG(pc,num_levels,NULL);CHKERRQ(ierr);
1718a2c336bSFande Kong   /* Set coarse matrices and interpolations to PCMG */
172360ee056SFande Kong   for (level=num_levels-1; level>0; level--) {
173360ee056SFande Kong     Mat P=0, pmat=0;
174360ee056SFande Kong     Vec b, x,r;
1758a2c336bSFande Kong     if (hmg->subcoarsening) {
1768a2c336bSFande Kong       /* Grow interpolation. In the future, we should use MAIJ */
1778a2c336bSFande Kong       ierr = PCHMGExpandInterpolation_Private(interpolations[level-1],&P,blocksize);CHKERRQ(ierr);
1788a2c336bSFande Kong       ierr = MatDestroy(&interpolations[level-1]);CHKERRQ(ierr);
179360ee056SFande Kong     } else {
1808a2c336bSFande Kong       P = interpolations[level-1];
181360ee056SFande Kong     }
182360ee056SFande Kong     ierr = MatCreateVecs(P,&b,&r);CHKERRQ(ierr);
183360ee056SFande Kong     ierr = PCMGSetInterpolation(pc,level,P);CHKERRQ(ierr);
184360ee056SFande Kong     ierr = PCMGSetRestriction(pc,level,P);CHKERRQ(ierr);
185360ee056SFande Kong     ierr = MatDestroy(&P);CHKERRQ(ierr);
1868a2c336bSFande Kong     /* We reuse the matrices when we do not do subspace coarsening */
1878a2c336bSFande Kong     if ((level-1)>=0 && !hmg->subcoarsening) {
1888a2c336bSFande Kong       pmat = operators[level-1];
1898a2c336bSFande Kong       ierr = PCMGSetOperators(pc,level-1,pmat,pmat);CHKERRQ(ierr);
190360ee056SFande Kong       ierr = MatDestroy(&pmat);CHKERRQ(ierr);
191360ee056SFande Kong     }
192360ee056SFande Kong     ierr = PCMGSetRhs(pc,level-1,b);CHKERRQ(ierr);
193360ee056SFande Kong 
194360ee056SFande Kong     ierr = PCMGSetR(pc,level,r);CHKERRQ(ierr);
195360ee056SFande Kong     ierr = VecDestroy(&r);CHKERRQ(ierr);
196360ee056SFande Kong 
197fd2dd295SFande Kong     ierr = VecDuplicate(b,&x);CHKERRQ(ierr);
198360ee056SFande Kong     ierr = PCMGSetX(pc,level-1,x);CHKERRQ(ierr);
199360ee056SFande Kong     ierr = VecDestroy(&x);CHKERRQ(ierr);
200360ee056SFande Kong     ierr = VecDestroy(&b);CHKERRQ(ierr);
201360ee056SFande Kong   }
2028a2c336bSFande Kong   ierr = PetscFree(interpolations);CHKERRQ(ierr);
2038a2c336bSFande Kong   if (!hmg->subcoarsening) {
2048a2c336bSFande Kong     ierr = PetscFree(operators);CHKERRQ(ierr);
2058a2c336bSFande Kong   }
2068a2c336bSFande Kong   /* Turn Galerkin off when we already have coarse operators */
2078a2c336bSFande Kong   ierr = PCMGSetGalerkin(pc,hmg->subcoarsening ? PC_MG_GALERKIN_PMAT:PC_MG_GALERKIN_NONE);CHKERRQ(ierr);
208360ee056SFande Kong   ierr = PCSetDM(pc,NULL);CHKERRQ(ierr);
209360ee056SFande Kong   ierr = PCSetUseAmat(pc,PETSC_FALSE);CHKERRQ(ierr);
210360ee056SFande Kong   ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
211360ee056SFande Kong   ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */
212360ee056SFande Kong   ierr = PetscOptionsEnd();CHKERRQ(ierr);
213360ee056SFande Kong   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
214360ee056SFande Kong   PetscFunctionReturn(0);
215360ee056SFande Kong }
216360ee056SFande Kong 
217360ee056SFande Kong PetscErrorCode PCDestroy_HMG(PC pc)
218360ee056SFande Kong {
219360ee056SFande Kong   PetscErrorCode ierr;
220360ee056SFande Kong   PC_MG          *mg  = (PC_MG*)pc->data;
2218a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
222360ee056SFande Kong 
223360ee056SFande Kong   PetscFunctionBegin;
2248a2c336bSFande Kong   ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr);
2258a2c336bSFande Kong   ierr = PetscFree(hmg->innerpctype);CHKERRQ(ierr);
2268a2c336bSFande Kong   ierr = PetscFree(hmg);CHKERRQ(ierr);
227360ee056SFande Kong   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
228fd2dd295SFande Kong 
229fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",NULL);CHKERRQ(ierr);
230fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",NULL);CHKERRQ(ierr);
231fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",NULL);CHKERRQ(ierr);
232360ee056SFande Kong   PetscFunctionReturn(0);
233360ee056SFande Kong }
234360ee056SFande Kong 
235360ee056SFande Kong PetscErrorCode PCView_HMG(PC pc,PetscViewer viewer)
236360ee056SFande Kong {
237360ee056SFande Kong   PC_MG          *mg = (PC_MG*)pc->data;
2388a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
239360ee056SFande Kong   PetscErrorCode ierr;
240360ee056SFande Kong   PetscBool      iascii;
241360ee056SFande Kong 
242360ee056SFande Kong   PetscFunctionBegin;
243360ee056SFande Kong   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
244360ee056SFande Kong   if (iascii) {
2458a2c336bSFande Kong     ierr = PetscViewerASCIIPrintf(viewer," Reuse interpolation: %s\n",hmg->reuseinterp? "true":"false");CHKERRQ(ierr);
2468a2c336bSFande Kong     ierr = PetscViewerASCIIPrintf(viewer," Use subspace coarsening: %s\n",hmg->subcoarsening? "true":"false");CHKERRQ(ierr);
2478a2c336bSFande Kong     ierr = PetscViewerASCIIPrintf(viewer," Inner PC type: %s \n",hmg->innerpctype);CHKERRQ(ierr);
248360ee056SFande Kong   }
249360ee056SFande Kong   ierr = PCView_MG(pc,viewer);CHKERRQ(ierr);
250360ee056SFande Kong   PetscFunctionReturn(0);
251360ee056SFande Kong }
252360ee056SFande Kong 
253360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems *PetscOptionsObject,PC pc)
254360ee056SFande Kong {
255360ee056SFande Kong   PetscErrorCode ierr;
256360ee056SFande Kong   PC_MG          *mg = (PC_MG*)pc->data;
2578a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
258360ee056SFande Kong 
259360ee056SFande Kong   PetscFunctionBegin;
260360ee056SFande Kong   ierr = PetscOptionsHead(PetscOptionsObject,"HMG");CHKERRQ(ierr);
2618a2c336bSFande 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);
2628a2c336bSFande Kong   ierr = PetscOptionsBool("-pc_hmg_use_subspace_coarsening","Use the subspace coarsening to compute the interpolations","None",hmg->subcoarsening,&hmg->subcoarsening,NULL);CHKERRQ(ierr);
263360ee056SFande Kong   ierr = PetscOptionsTail();CHKERRQ(ierr);
264360ee056SFande Kong   PetscFunctionReturn(0);
265360ee056SFande Kong }
266360ee056SFande Kong 
267fd2dd295SFande Kong static PetscErrorCode PCHMGSetReuseInterpolation_HMG(PC pc, PetscBool reuse)
268fd2dd295SFande Kong {
269*07a4832bSFande Kong   PC_MG          *mg  = (PC_MG*)pc->data;
270*07a4832bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
271fd2dd295SFande Kong 
272fd2dd295SFande Kong   PetscFunctionBegin;
273fd2dd295SFande Kong   hmg->reuseinterp = reuse;
274fd2dd295SFande Kong   PetscFunctionReturn(0);
275fd2dd295SFande Kong }
276fd2dd295SFande Kong 
277360ee056SFande Kong /*MC
2788a2c336bSFande Kong    PCHMGSetReuseInterpolation - Reuse interpolation matrices in HMG
2798a2c336bSFande Kong 
2808a2c336bSFande Kong    Logically Collective on PC
2818a2c336bSFande Kong 
2828a2c336bSFande Kong    Input Parameters:
2838a2c336bSFande Kong +  pc - the HMG context
2848a2c336bSFande Kong -  reuse - True indicates that HMG will reuse the interpolations
2858a2c336bSFande Kong 
2868a2c336bSFande Kong    Options Database Keys:
2878a2c336bSFande Kong +  -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
2888a2c336bSFande Kong 
2898a2c336bSFande Kong    Level: beginner
2908a2c336bSFande Kong 
2918a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, reuse, set
2928a2c336bSFande Kong 
2938a2c336bSFande Kong .seealso: PCHMG
2948a2c336bSFande Kong M*/
2958a2c336bSFande Kong PetscErrorCode PCHMGSetReuseInterpolation(PC pc, PetscBool reuse)
2968a2c336bSFande Kong {
297fd2dd295SFande Kong   PetscErrorCode ierr;
298fd2dd295SFande Kong 
299fd2dd295SFande Kong   PetscFunctionBegin;
300fd2dd295SFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
301fd2dd295SFande Kong   ierr = PetscUseMethod(pc,"PCHMGSetReuseInterpolation_C",(PC,PetscBool),(pc,reuse));CHKERRQ(ierr);
302fd2dd295SFande Kong   PetscFunctionReturn(0);
303fd2dd295SFande Kong }
304fd2dd295SFande Kong 
305fd2dd295SFande Kong static PetscErrorCode PCHMGSetUseSubspaceCoarsening_HMG(PC pc, PetscBool subspace)
306fd2dd295SFande Kong {
307*07a4832bSFande Kong   PC_MG          *mg  = (PC_MG*)pc->data;
308*07a4832bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
3098a2c336bSFande Kong 
3108a2c336bSFande Kong   PetscFunctionBegin;
311fd2dd295SFande Kong   hmg->subcoarsening = subspace;
3128a2c336bSFande Kong   PetscFunctionReturn(0);
3138a2c336bSFande Kong }
3148a2c336bSFande Kong 
3158a2c336bSFande Kong /*MC
3168a2c336bSFande Kong    PCHMGSetUseSubspaceCoarsening - Use subspace coarsening in HMG
3178a2c336bSFande Kong 
3188a2c336bSFande Kong    Logically Collective on PC
3198a2c336bSFande Kong 
3208a2c336bSFande Kong    Input Parameters:
3218a2c336bSFande Kong +  pc - the HMG context
3228a2c336bSFande Kong -  reuse - True indicates that HMG will use the subspace coarsening
3238a2c336bSFande Kong 
3248a2c336bSFande Kong    Options Database Keys:
3258a2c336bSFande Kong +  -pc_hmg_use_subspace_coarsening  <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix).
3268a2c336bSFande Kong 
3278a2c336bSFande Kong    Level: beginner
3288a2c336bSFande Kong 
3298a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, subspace, coarsening
3308a2c336bSFande Kong 
3318a2c336bSFande Kong .seealso: PCHMG
3328a2c336bSFande Kong M*/
3338a2c336bSFande Kong PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC pc, PetscBool subspace)
3348a2c336bSFande Kong {
335fd2dd295SFande Kong   PetscErrorCode ierr;
3368a2c336bSFande Kong 
3378a2c336bSFande Kong   PetscFunctionBegin;
3388a2c336bSFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
339fd2dd295SFande Kong   ierr = PetscUseMethod(pc,"PCHMGSetUseSubspaceCoarsening_C",(PC,PetscBool),(pc,subspace));CHKERRQ(ierr);
340fd2dd295SFande Kong   PetscFunctionReturn(0);
341fd2dd295SFande Kong }
342fd2dd295SFande Kong 
343fd2dd295SFande Kong static PetscErrorCode PCHMGSetInnerPCType_HMG(PC pc, PCType type)
344fd2dd295SFande Kong {
345*07a4832bSFande Kong   PC_MG           *mg  = (PC_MG*)pc->data;
346*07a4832bSFande Kong   PC_HMG          *hmg = (PC_HMG*) mg->innerctx;
347fd2dd295SFande Kong   PetscErrorCode  ierr;
348fd2dd295SFande Kong 
349fd2dd295SFande Kong   PetscFunctionBegin;
350fd2dd295SFande Kong   ierr = PetscStrallocpy(type,&(hmg->innerpctype));CHKERRQ(ierr);
3518a2c336bSFande Kong   PetscFunctionReturn(0);
3528a2c336bSFande Kong }
3538a2c336bSFande Kong 
3548a2c336bSFande Kong /*MC
3558a2c336bSFande Kong    PCHMGSetInnerPCType - Set an inner PC type
3568a2c336bSFande Kong 
3578a2c336bSFande Kong    Logically Collective on PC
3588a2c336bSFande Kong 
3598a2c336bSFande Kong    Input Parameters:
3608a2c336bSFande Kong +  pc - the HMG context
3618a2c336bSFande Kong -  type - <hypre, gamg> coarsening algorithm
3628a2c336bSFande Kong 
3638a2c336bSFande Kong    Options Database Keys:
364fd2dd295SFande Kong +  -hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix
3658a2c336bSFande Kong 
3668a2c336bSFande Kong    Level: beginner
3678a2c336bSFande Kong 
3688a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, coarsening
3698a2c336bSFande Kong 
3708a2c336bSFande Kong .seealso: PCHMG, PCType
3718a2c336bSFande Kong M*/
3728a2c336bSFande Kong PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type)
3738a2c336bSFande Kong {
3748a2c336bSFande Kong   PetscErrorCode  ierr;
3758a2c336bSFande Kong 
3768a2c336bSFande Kong   PetscFunctionBegin;
3778a2c336bSFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
378fd2dd295SFande Kong   ierr = PetscUseMethod(pc,"PCHMGSetInnerPCType_C",(PC,PCType),(pc,type));CHKERRQ(ierr);
3798a2c336bSFande Kong   PetscFunctionReturn(0);
3808a2c336bSFande Kong }
3818a2c336bSFande Kong 
3828a2c336bSFande Kong /*MC
383fd2dd295SFande Kong    PCHMG - Hybrid of PETSc preconditioners (such as ASM, BJacobi, SOR, etc.) and Hypre BoomerAMG, GAMG or other multilevel methods. BoomerAMG, GAMG
384fd2dd295SFande Kong            or other multilevel methods is used to coarsen matrix and generate a sequence of coarse matrices and interpolations. The matrices and
385fd2dd295SFande Kong            interpolations are employed to construct PCMG, and then any available PETSc preconditioners can be chosen as smoothers and the coarse solver.
386360ee056SFande Kong 
387360ee056SFande Kong    Options Database Keys:
3888a2c336bSFande Kong +  -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
3898a2c336bSFande Kong .  -pc_hmg_use_subspace_coarsening  <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix).
390fd2dd295SFande Kong .  -hmg_inner_pc_type <hypre, gamg, ...> - What method is used to coarsen matrix
391360ee056SFande Kong 
392360ee056SFande Kong 
393360ee056SFande Kong    Notes:
394360ee056SFande Kong     For multicomponent problems, we can just coarsen one submatrix associated with one particular component. In this way, the preconditioner setup
395360ee056SFande Kong     time is significantly reduced. One typical use case is neutron transport equations. There are many variables on each mesh vertex due to the
396360ee056SFande Kong     of angle and energy. Each variable, in fact, corresponds to the same PDEs but with different material properties.
397360ee056SFande Kong 
398360ee056SFande Kong    Level: beginner
399360ee056SFande Kong 
4008a2c336bSFande Kong    Concepts: Hybrid of ASM and MG, Subspace Coarsening
401360ee056SFande Kong 
402360ee056SFande Kong     References:
403360ee056SFande 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
404360ee056SFande Kong     Newton-Krylov-Schwarz method with subspace-based coarsening and partition-based balancing for the multigroup neutron transport equations on
405360ee056SFande Kong     3D unstructured meshes, arXiv preprint arXiv:1903.03659, 2019
406360ee056SFande Kong 
407fd2dd295SFande Kong .seealso:  PCCreate(), PCSetType(), PCType, PC, PCMG, PCHYPRE, PCHMG, PCGetCoarseOperators(), PCGetInterpolations(), PCHMGSetReuseInterpolation(), PCHMGSetUseSubspaceCoarsening(),
408fd2dd295SFande Kong            PCHMGSetInnerPCType()
409360ee056SFande Kong 
410360ee056SFande Kong M*/
411360ee056SFande Kong PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc)
412360ee056SFande Kong {
413360ee056SFande Kong   PetscErrorCode ierr;
414360ee056SFande Kong   PC_HMG         *hmg;
415360ee056SFande Kong   PC_MG          *mg;
416360ee056SFande Kong 
417360ee056SFande Kong   PetscFunctionBegin;
418360ee056SFande Kong   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
419360ee056SFande Kong   if (pc->ops->destroy) {
420360ee056SFande Kong     ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr);
421360ee056SFande Kong     pc->data = 0;
422360ee056SFande Kong   }
423360ee056SFande Kong   ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr);
424360ee056SFande Kong 
425360ee056SFande Kong   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr);
426fd2dd295SFande Kong   ierr = PetscNew(&hmg);CHKERRQ(ierr);
427360ee056SFande Kong 
428360ee056SFande Kong   mg                      = (PC_MG*) pc->data;
429360ee056SFande Kong   mg->innerctx            = hmg;
430360ee056SFande Kong   hmg->reuseinterp        = PETSC_FALSE;
4318a2c336bSFande Kong   hmg->subcoarsening      = PETSC_FALSE;
4328a2c336bSFande Kong   hmg->innerpc            = NULL;
433360ee056SFande Kong 
434360ee056SFande Kong   pc->ops->setfromoptions = PCSetFromOptions_HMG;
435360ee056SFande Kong   pc->ops->view           = PCView_HMG;
436360ee056SFande Kong   pc->ops->destroy        = PCDestroy_HMG;
437360ee056SFande Kong   pc->ops->setup          = PCSetUp_HMG;
438fd2dd295SFande Kong 
439fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",PCHMGSetReuseInterpolation_HMG);CHKERRQ(ierr);
440fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",PCHMGSetUseSubspaceCoarsening_HMG);CHKERRQ(ierr);
441fd2dd295SFande Kong   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",PCHMGSetInnerPCType_HMG);CHKERRQ(ierr);
442360ee056SFande Kong   PetscFunctionReturn(0);
443360ee056SFande Kong }
444