xref: /petsc/src/ksp/pc/impls/hmg/hmg.c (revision 8a2c336b4146966a0f1271f829a3a6376f16d9c3)
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 {
8*8a2c336bSFande Kong   PC               innerpc;            /* A MG inner PC (Hypre or PCGAMG) to setup interpolations and coarse operators  */
9*8a2c336bSFande Kong   char*            innerpctype;        /* PCGAMG or PCHYPRE */
10*8a2c336bSFande Kong   PetscBool        reuseinterp;        /* A flag indicates if or not to reuse the interpolations */
11*8a2c336bSFande Kong   PetscBool        subcoarsening;
12*8a2c336bSFande Kong   PetscErrorCode   (*getoperators)(PC,PetscInt*,Mat**);
13*8a2c336bSFande Kong   PetscErrorCode   (*getinterpolations)(PC,PetscInt*,Mat**);
14*8a2c336bSFande Kong   PetscErrorCode   (*setfromoptions)(PetscOptionItems*,PC);
15360ee056SFande Kong } PC_HMG;
16360ee056SFande Kong 
17360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems*,PC);
18360ee056SFande Kong PetscErrorCode PCReset_MG(PC);
19*8a2c336bSFande Kong #if PETSC_HAVE_HYPRE
20*8a2c336bSFande Kong PetscErrorCode PCHYPREBoomerAMGGetCoarseOperators(PC,PetscInt*,Mat**);
21*8a2c336bSFande Kong PetscErrorCode PCHYPREBoomerAMGGetInterpolations(PC,PetscInt*,Mat**);
22*8a2c336bSFande Kong PetscErrorCode PCSetFromOptions_HYPRE(PetscOptionItems*,PC);
23*8a2c336bSFande Kong #endif
24*8a2c336bSFande Kong PetscErrorCode PCMGGetCoarseOperators(PC,PetscInt*,Mat**);
25*8a2c336bSFande Kong PetscErrorCode PCMGGetInterpolations(PC,PetscInt*,Mat **);
26*8a2c336bSFande Kong PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems*,PC);
27360ee056SFande Kong 
28*8a2c336bSFande Kong static PetscErrorCode PCHMGSetInnerPCFromOptions_Private(PC pc,PetscOptionItems *PetscOptionsObject)
29*8a2c336bSFande Kong {
30*8a2c336bSFande Kong   PetscErrorCode     ierr;
31*8a2c336bSFande Kong   PC_MG              *mg   = (PC_MG*)pc->data;
32*8a2c336bSFande Kong   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
33*8a2c336bSFande Kong 
34*8a2c336bSFande Kong   PetscFunctionBegin;
35*8a2c336bSFande Kong   if (hmg->setfromoptions) {
36*8a2c336bSFande Kong     ierr = (*hmg->setfromoptions)(PetscOptionsObject,hmg->innerpc);CHKERRQ(ierr);
37*8a2c336bSFande Kong   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_NULL,"Inner PC does not support SetFromOptions\n");
38*8a2c336bSFande Kong   PetscFunctionReturn(0);
39*8a2c336bSFande Kong }
40*8a2c336bSFande Kong 
41*8a2c336bSFande Kong static PetscErrorCode PCHMGGetInnerPCCoarseOperators_Private(PC pc,PetscInt *nlevels,Mat **operators)
42*8a2c336bSFande Kong {
43*8a2c336bSFande Kong   PetscErrorCode     ierr;
44*8a2c336bSFande Kong   PC_MG              *mg   = (PC_MG*)pc->data;
45*8a2c336bSFande Kong   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
46*8a2c336bSFande Kong 
47*8a2c336bSFande Kong   PetscFunctionBegin;
48*8a2c336bSFande Kong   if (hmg->getoperators) {
49*8a2c336bSFande Kong     ierr = (*hmg->getoperators)(hmg->innerpc,nlevels,operators);CHKERRQ(ierr);
50*8a2c336bSFande Kong   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_NULL,"Inner PC does not support GetOperators\n");
51*8a2c336bSFande Kong   PetscFunctionReturn(0);
52*8a2c336bSFande Kong }
53*8a2c336bSFande Kong 
54*8a2c336bSFande Kong static PetscErrorCode PCHMGGetInnerPCInterpolations_Private(PC pc,PetscInt *nlevels,Mat *interpolations[])
55*8a2c336bSFande Kong {
56*8a2c336bSFande Kong   PetscErrorCode     ierr;
57*8a2c336bSFande Kong   PC_MG              *mg   = (PC_MG*)pc->data;
58*8a2c336bSFande Kong   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
59*8a2c336bSFande Kong 
60*8a2c336bSFande Kong   PetscFunctionBegin;
61*8a2c336bSFande Kong   if (hmg->getinterpolations) {
62*8a2c336bSFande Kong     ierr = (*hmg->getinterpolations)(hmg->innerpc,nlevels,interpolations);CHKERRQ(ierr);
63*8a2c336bSFande Kong   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_NULL,"Inner PC does not support GetInterpolation\n");
64*8a2c336bSFande Kong   PetscFunctionReturn(0);
65*8a2c336bSFande Kong }
66*8a2c336bSFande Kong 
67*8a2c336bSFande Kong static PetscErrorCode PCHMGExtractSubMatrix_Private(Mat pmat,Mat *submat,MatReuse reuse,PetscInt blocksize)
68360ee056SFande Kong {
69360ee056SFande Kong   IS             isrow;
70360ee056SFande Kong   PetscErrorCode ierr;
71*8a2c336bSFande Kong   PetscInt       rstart,rend;
72360ee056SFande Kong   MPI_Comm       comm;
73360ee056SFande Kong 
74360ee056SFande Kong   PetscFunctionBegin;
75360ee056SFande Kong   ierr = PetscObjectGetComm((PetscObject)pmat,&comm);CHKERRQ(ierr);
76360ee056SFande Kong   ierr = MatGetOwnershipRange(pmat,&rstart,&rend);CHKERRQ(ierr);
77360ee056SFande Kong   if ((rend-rstart)%blocksize != 0) SETERRQ3(comm,PETSC_ERR_ARG_INCOMP,"Block size %d is inconsisent for [%d, %d) \n",blocksize,rstart,rend);
78*8a2c336bSFande Kong   ierr = ISCreateStride(comm,(rend-rstart)/blocksize,rstart,blocksize,&isrow);
79360ee056SFande Kong   ierr = MatCreateSubMatrix(pmat,isrow,isrow,reuse,submat);CHKERRQ(ierr);
80360ee056SFande Kong   ierr = ISDestroy(&isrow);CHKERRQ(ierr);
81360ee056SFande Kong   PetscFunctionReturn(0);
82360ee056SFande Kong }
83360ee056SFande Kong 
84*8a2c336bSFande Kong static PetscErrorCode PCHMGExpandInterpolation_Private(Mat subinterp, Mat *interp, PetscInt blocksize)
85360ee056SFande Kong {
86360ee056SFande Kong   PetscInt              subrstart,subrend,subrowsize,subcolsize,subcstart,subcend,rowsize,colsize;
87*8a2c336bSFande Kong   PetscInt              subrow,row,nz,*d_nnz,*o_nnz,i,j,dnz,onz,max_nz,*indices;
88*8a2c336bSFande Kong   const PetscInt        *idx;
89*8a2c336bSFande Kong   const PetscScalar     *values;
90360ee056SFande Kong   PetscErrorCode        ierr;
91*8a2c336bSFande Kong   MPI_Comm              comm;
92360ee056SFande Kong 
93360ee056SFande Kong   PetscFunctionBegin;
94*8a2c336bSFande Kong   ierr = PetscObjectGetComm((PetscObject)subinterp,&comm);CHKERRQ(ierr);
95*8a2c336bSFande Kong   ierr = MatGetOwnershipRange(subinterp,&subrstart,&subrend);CHKERRQ(ierr);
96*8a2c336bSFande Kong   subrowsize = subrend-subrstart;
97360ee056SFande Kong   rowsize = subrowsize*blocksize;
98360ee056SFande Kong   ierr = PetscCalloc2(rowsize,&d_nnz,rowsize,&o_nnz);CHKERRQ(ierr);
99*8a2c336bSFande Kong   ierr = MatGetOwnershipRangeColumn(subinterp,&subcstart,&subcend);CHKERRQ(ierr);
100*8a2c336bSFande Kong   subcolsize = subcend - subcstart;
101*8a2c336bSFande Kong   colsize    = subcolsize*blocksize;
102360ee056SFande Kong   max_nz = 0;
103*8a2c336bSFande Kong   for (subrow=subrstart;subrow<subrend;subrow++) {
104*8a2c336bSFande Kong     ierr = MatGetRow(subinterp,subrow,&nz,&idx,NULL);CHKERRQ(ierr);
105360ee056SFande Kong     if (max_nz<nz) max_nz = nz;
106360ee056SFande Kong     dnz = 0; onz = 0;
107360ee056SFande Kong     for (i=0;i<nz;i++) {
108*8a2c336bSFande Kong       if(idx[i]>=subcstart && idx[i]<subcend) dnz++;
109*8a2c336bSFande Kong       else onz++;
110360ee056SFande Kong     }
111360ee056SFande Kong     for (i=0;i<blocksize;i++) {
112360ee056SFande Kong       d_nnz[(subrow-subrstart)*blocksize+i] = dnz;
113360ee056SFande Kong       o_nnz[(subrow-subrstart)*blocksize+i] = onz;
114360ee056SFande Kong     }
115*8a2c336bSFande Kong     ierr = MatRestoreRow(subinterp,subrow,&nz,&idx,NULL);CHKERRQ(ierr);
116360ee056SFande Kong   }
117360ee056SFande Kong   ierr = MatCreateAIJ(comm,rowsize,colsize,PETSC_DETERMINE,PETSC_DETERMINE,0,d_nnz,0,o_nnz,interp);CHKERRQ(ierr);
118360ee056SFande Kong   ierr = MatSetOption(*interp,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
119360ee056SFande Kong   ierr = MatSetOption(*interp,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
120360ee056SFande Kong   ierr = MatSetOption(*interp,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
121360ee056SFande Kong   ierr = MatSetFromOptions(*interp);CHKERRQ(ierr);
122360ee056SFande Kong 
123360ee056SFande Kong   ierr = MatSetUp(*interp);CHKERRQ(ierr);
124360ee056SFande Kong   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
125360ee056SFande Kong   ierr = PetscCalloc1(max_nz,&indices);CHKERRQ(ierr);
126*8a2c336bSFande Kong   for (subrow=subrstart; subrow<subrend; subrow++) {
127*8a2c336bSFande Kong     ierr = MatGetRow(subinterp,subrow,&nz,&idx,&values);CHKERRQ(ierr);
128360ee056SFande Kong     for (i=0;i<blocksize;i++) {
129360ee056SFande Kong       row = subrow*blocksize+i;
130360ee056SFande Kong       for (j=0;j<nz;j++) {
131360ee056SFande Kong         indices[j] = idx[j]*blocksize+i;
132360ee056SFande Kong       }
133360ee056SFande Kong       ierr = MatSetValues(*interp,1,&row,nz,indices,values,INSERT_VALUES);CHKERRQ(ierr);
134360ee056SFande Kong     }
135*8a2c336bSFande Kong     ierr = MatRestoreRow(subinterp,subrow,&nz,&idx,&values);CHKERRQ(ierr);
136360ee056SFande Kong   }
137360ee056SFande Kong   ierr = PetscFree(indices);CHKERRQ(ierr);
138360ee056SFande Kong   ierr = MatAssemblyBegin(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139360ee056SFande Kong   ierr = MatAssemblyEnd(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
140360ee056SFande Kong   PetscFunctionReturn(0);
141360ee056SFande Kong }
142360ee056SFande Kong 
143360ee056SFande Kong PetscErrorCode PCSetUp_HMG(PC pc)
144360ee056SFande Kong {
145360ee056SFande Kong   PetscErrorCode     ierr;
146360ee056SFande Kong   Mat                PA, submat;
147360ee056SFande Kong   PC_MG              *mg   = (PC_MG*)pc->data;
148360ee056SFande Kong   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
149360ee056SFande Kong   MPI_Comm           comm;
150360ee056SFande Kong   PetscInt           level;
151360ee056SFande Kong   PetscInt           num_levels;
152*8a2c336bSFande Kong   PetscBool          same;
153*8a2c336bSFande Kong   Mat                *operators,*interpolations;
154*8a2c336bSFande Kong   PetscInt           blocksize;
155360ee056SFande Kong 
156360ee056SFande Kong   PetscFunctionBegin;
157360ee056SFande Kong   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
158360ee056SFande Kong   if (pc->setupcalled) {
159*8a2c336bSFande Kong    /* Only reuse interpolations when nonzero pattern does not change */
160360ee056SFande Kong    if (pc->flag == SAME_NONZERO_PATTERN && hmg->reuseinterp) {
161360ee056SFande Kong      ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT);CHKERRQ(ierr);
162360ee056SFande Kong      ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
163360ee056SFande Kong      PetscFunctionReturn(0);
164360ee056SFande Kong     } else {
165360ee056SFande Kong      ierr = PCReset_MG(pc);CHKERRQ(ierr);
166360ee056SFande Kong      pc->setupcalled = PETSC_FALSE;
167360ee056SFande Kong     }
168360ee056SFande Kong   }
169360ee056SFande Kong 
170*8a2c336bSFande Kong   /* Create an inner PC (GAMG or HYPRE) */
171*8a2c336bSFande Kong   if (!hmg->innerpc) {
172*8a2c336bSFande Kong     ierr = PCCreate(comm,&hmg->innerpc);CHKERRQ(ierr);
173*8a2c336bSFande Kong     ierr = PCSetType(hmg->innerpc,hmg->innerpctype);CHKERRQ(ierr);
174*8a2c336bSFande Kong     ierr = PetscStrcmp(hmg->innerpctype,"hypre",&same);CHKERRQ(ierr);
175*8a2c336bSFande Kong     if (same) {
176*8a2c336bSFande Kong       ierr = PCHYPRESetType(hmg->innerpc,"boomeramg");CHKERRQ(ierr);
177*8a2c336bSFande Kong     }
178360ee056SFande Kong   }
179360ee056SFande Kong   ierr = PCGetOperators(pc,NULL,&PA);CHKERRQ(ierr);
180*8a2c336bSFande Kong   /* Users need to correctly set a block size of matrix in order to use subspace coarsening */
181*8a2c336bSFande Kong   ierr = MatGetBlockSize(PA,&blocksize);CHKERRQ(ierr);
182*8a2c336bSFande Kong   if (blocksize<=1) hmg->subcoarsening = PETSC_FALSE;
183*8a2c336bSFande Kong   /* Extract a submatrix for constructing subinterpolations */
184*8a2c336bSFande Kong   if (hmg->subcoarsening) {
185*8a2c336bSFande Kong     ierr = PCHMGExtractSubMatrix_Private(PA,&submat,MAT_INITIAL_MATRIX,blocksize);CHKERRQ(ierr);
186360ee056SFande Kong     PA = submat;
187360ee056SFande Kong   }
188*8a2c336bSFande Kong   ierr = PCSetOperators(hmg->innerpc,PA,PA);CHKERRQ(ierr);
189*8a2c336bSFande Kong   if (hmg->subcoarsening) {
190360ee056SFande Kong    ierr = MatDestroy(&PA);CHKERRQ(ierr);
191360ee056SFande Kong   }
192*8a2c336bSFande Kong   /* Setup inner PC correctly. During this step, matrix will be coarsened */
193*8a2c336bSFande Kong   ierr = PCSetUseAmat(hmg->innerpc,PETSC_FALSE);CHKERRQ(ierr);
194*8a2c336bSFande Kong   ierr = PetscObjectOptionsBegin((PetscObject)hmg->innerpc);CHKERRQ(ierr);
195*8a2c336bSFande Kong   ierr = PCHMGSetInnerPCFromOptions_Private(pc,PetscOptionsObject);CHKERRQ(ierr);
196360ee056SFande Kong   ierr = PetscOptionsEnd();CHKERRQ(ierr);
197*8a2c336bSFande Kong   ierr = PCSetUp(hmg->innerpc);
198360ee056SFande Kong 
199*8a2c336bSFande Kong   /* Obtain interpolations IN PLACE. For BoomerAMG, (I,J,data) is reused to avoid memory overhead */
200*8a2c336bSFande Kong   ierr = PCHMGGetInnerPCInterpolations_Private(pc,&num_levels,&interpolations);CHKERRQ(ierr);
201*8a2c336bSFande Kong   /* We can reuse the coarse operators when we do the full space coarsening */
202*8a2c336bSFande Kong   if (!hmg->subcoarsening) {
203*8a2c336bSFande Kong     ierr = PCHMGGetInnerPCCoarseOperators_Private(pc,&num_levels,&operators);CHKERRQ(ierr);
204360ee056SFande Kong   }
205360ee056SFande Kong 
206*8a2c336bSFande Kong   ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr);
207*8a2c336bSFande Kong   hmg->innerpc = 0;
208*8a2c336bSFande Kong   ierr = PCMGSetLevels_MG(pc,num_levels,NULL);CHKERRQ(ierr);
209*8a2c336bSFande Kong   /* Set coarse matrices and interpolations to PCMG */
210360ee056SFande Kong   for (level=num_levels-1; level>0; level--) {
211360ee056SFande Kong     Mat P=0, pmat=0;
212360ee056SFande Kong     Vec b, x,r;
213*8a2c336bSFande Kong     if (hmg->subcoarsening) {
214*8a2c336bSFande Kong       /* Grow interpolation. In the future, we should use MAIJ */
215*8a2c336bSFande Kong       ierr = PCHMGExpandInterpolation_Private(interpolations[level-1],&P,blocksize);CHKERRQ(ierr);
216*8a2c336bSFande Kong       ierr = MatDestroy(&interpolations[level-1]);CHKERRQ(ierr);
217360ee056SFande Kong     } else {
218*8a2c336bSFande Kong       P = interpolations[level-1];
219360ee056SFande Kong     }
220360ee056SFande Kong     ierr = MatCreateVecs(P,&b,&r);CHKERRQ(ierr);
221360ee056SFande Kong     ierr = PCMGSetInterpolation(pc,level,P);CHKERRQ(ierr);
222360ee056SFande Kong     ierr = PCMGSetRestriction(pc,level,P);CHKERRQ(ierr);
223360ee056SFande Kong     ierr = MatDestroy(&P);CHKERRQ(ierr);
224*8a2c336bSFande Kong     /* We reuse the matrices when we do not do subspace coarsening */
225*8a2c336bSFande Kong     if ((level-1)>=0 && !hmg->subcoarsening) {
226*8a2c336bSFande Kong       pmat = operators[level-1];
227*8a2c336bSFande Kong       ierr = PCMGSetOperators(pc,level-1,pmat,pmat);CHKERRQ(ierr);
228360ee056SFande Kong       ierr = MatDestroy(&pmat);CHKERRQ(ierr);
229360ee056SFande Kong     }
230360ee056SFande Kong     ierr = PCMGSetRhs(pc,level-1,b);CHKERRQ(ierr);
231360ee056SFande Kong 
232360ee056SFande Kong     ierr = PCMGSetR(pc,level,r);CHKERRQ(ierr);
233360ee056SFande Kong     ierr = VecDestroy(&r);CHKERRQ(ierr);
234360ee056SFande Kong 
235360ee056SFande Kong     ierr = VecDuplicate(b,&x);
236360ee056SFande Kong     ierr = PCMGSetX(pc,level-1,x);CHKERRQ(ierr);
237360ee056SFande Kong     ierr = VecDestroy(&x);CHKERRQ(ierr);
238360ee056SFande Kong     ierr = VecDestroy(&b);CHKERRQ(ierr);
239360ee056SFande Kong   }
240*8a2c336bSFande Kong   ierr = PetscFree(interpolations);CHKERRQ(ierr);
241*8a2c336bSFande Kong   if (!hmg->subcoarsening) {
242*8a2c336bSFande Kong     ierr = PetscFree(operators);CHKERRQ(ierr);
243*8a2c336bSFande Kong   }
244*8a2c336bSFande Kong   /* Turn Galerkin off when we already have coarse operators */
245*8a2c336bSFande Kong   ierr = PCMGSetGalerkin(pc,hmg->subcoarsening ? PC_MG_GALERKIN_PMAT:PC_MG_GALERKIN_NONE);CHKERRQ(ierr);
246360ee056SFande Kong   ierr = PCSetDM(pc,NULL);CHKERRQ(ierr);
247360ee056SFande Kong   ierr = PCSetUseAmat(pc,PETSC_FALSE);CHKERRQ(ierr);
248360ee056SFande Kong   ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
249360ee056SFande Kong   ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */
250360ee056SFande Kong   ierr = PetscOptionsEnd();CHKERRQ(ierr);
251360ee056SFande Kong   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
252360ee056SFande Kong   PetscFunctionReturn(0);
253360ee056SFande Kong }
254360ee056SFande Kong 
255360ee056SFande Kong PetscErrorCode PCDestroy_HMG(PC pc)
256360ee056SFande Kong {
257360ee056SFande Kong   PetscErrorCode ierr;
258360ee056SFande Kong   PC_MG          *mg  = (PC_MG*)pc->data;
259*8a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
260360ee056SFande Kong 
261360ee056SFande Kong   PetscFunctionBegin;
262*8a2c336bSFande Kong   ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr);
263*8a2c336bSFande Kong   ierr = PetscFree(hmg->innerpctype);CHKERRQ(ierr);
264*8a2c336bSFande Kong   ierr = PetscFree(hmg);CHKERRQ(ierr);
265360ee056SFande Kong   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
266360ee056SFande Kong   PetscFunctionReturn(0);
267360ee056SFande Kong }
268360ee056SFande Kong 
269360ee056SFande Kong PetscErrorCode PCView_HMG(PC pc,PetscViewer viewer)
270360ee056SFande Kong {
271360ee056SFande Kong   PC_MG          *mg = (PC_MG*)pc->data;
272*8a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
273360ee056SFande Kong   PetscErrorCode ierr;
274360ee056SFande Kong   PetscBool      iascii;
275360ee056SFande Kong 
276360ee056SFande Kong   PetscFunctionBegin;
277360ee056SFande Kong   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
278360ee056SFande Kong   if (iascii) {
279*8a2c336bSFande Kong     ierr = PetscViewerASCIIPrintf(viewer," Reuse interpolation: %s\n",hmg->reuseinterp? "true":"false");CHKERRQ(ierr);
280*8a2c336bSFande Kong     ierr = PetscViewerASCIIPrintf(viewer," Use subspace coarsening: %s\n",hmg->subcoarsening? "true":"false");CHKERRQ(ierr);
281*8a2c336bSFande Kong     ierr = PetscViewerASCIIPrintf(viewer," Inner PC type: %s \n",hmg->innerpctype);CHKERRQ(ierr);
282360ee056SFande Kong   }
283360ee056SFande Kong   ierr = PCView_MG(pc,viewer);CHKERRQ(ierr);
284360ee056SFande Kong   PetscFunctionReturn(0);
285360ee056SFande Kong }
286360ee056SFande Kong 
287*8a2c336bSFande Kong #if PETSC_HAVE_HYPRE
288*8a2c336bSFande Kong static const char *InnerPCType[]   = {"gamg","hypre"};
289*8a2c336bSFande Kong #else
290*8a2c336bSFande Kong static const char *InnerPCType[]   = {"gamg"};
291*8a2c336bSFande Kong #endif
292*8a2c336bSFande Kong 
293360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems *PetscOptionsObject,PC pc)
294360ee056SFande Kong {
295360ee056SFande Kong   PetscErrorCode ierr;
296360ee056SFande Kong   PC_MG          *mg = (PC_MG*)pc->data;
297*8a2c336bSFande Kong   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
298*8a2c336bSFande Kong   PetscBool      flg;
299*8a2c336bSFande Kong   PetscInt       indx;
300360ee056SFande Kong 
301360ee056SFande Kong   PetscFunctionBegin;
302360ee056SFande Kong   ierr = PetscOptionsHead(PetscOptionsObject,"HMG");CHKERRQ(ierr);
303*8a2c336bSFande 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);
304*8a2c336bSFande Kong   ierr = PetscOptionsBool("-pc_hmg_use_subspace_coarsening","Use the subspace coarsening to compute the interpolations","None",hmg->subcoarsening,&hmg->subcoarsening,NULL);CHKERRQ(ierr);
305*8a2c336bSFande Kong #if PETSC_HAVE_HYPRE
306*8a2c336bSFande Kong   indx = 1;
307*8a2c336bSFande Kong #else
308*8a2c336bSFande Kong   indx = 0;
309*8a2c336bSFande Kong #endif
310*8a2c336bSFande Kong   ierr = PetscOptionsEList("-pc_hmg_inner_pc_type","Inner PC type",NULL,InnerPCType,2,"gamg",&indx,&flg);CHKERRQ(ierr);
311*8a2c336bSFande Kong  if (indx == 0) {
312*8a2c336bSFande Kong    hmg->getinterpolations = PCMGGetInterpolations;
313*8a2c336bSFande Kong    hmg->getoperators      = PCMGGetCoarseOperators;
314*8a2c336bSFande Kong    hmg->setfromoptions    = PCSetFromOptions_GAMG;
315*8a2c336bSFande Kong    ierr = PetscStrallocpy(PCGAMG,&(hmg->innerpctype));CHKERRQ(ierr);
316*8a2c336bSFande Kong   }
317*8a2c336bSFande Kong #if PETSC_HAVE_HYPRE
318*8a2c336bSFande Kong  else if (indx == 1) {
319*8a2c336bSFande Kong     hmg->getinterpolations = PCHYPREBoomerAMGGetInterpolations;
320*8a2c336bSFande Kong     hmg->getoperators      = PCHYPREBoomerAMGGetCoarseOperators;
321*8a2c336bSFande Kong     hmg->setfromoptions    = PCSetFromOptions_HYPRE;
322*8a2c336bSFande Kong     ierr = PetscStrallocpy(PCHYPRE,&(hmg->innerpctype));CHKERRQ(ierr);
323*8a2c336bSFande Kong   }
324*8a2c336bSFande Kong #endif
325*8a2c336bSFande Kong   else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown inner PC type");
326*8a2c336bSFande Kong 
327360ee056SFande Kong   ierr = PetscOptionsTail();CHKERRQ(ierr);
328360ee056SFande Kong   PetscFunctionReturn(0);
329360ee056SFande Kong }
330360ee056SFande Kong 
331360ee056SFande Kong /*MC
332*8a2c336bSFande Kong    PCHMGSetReuseInterpolation - Reuse interpolation matrices in HMG
333*8a2c336bSFande Kong 
334*8a2c336bSFande Kong    Logically Collective on PC
335*8a2c336bSFande Kong 
336*8a2c336bSFande Kong    Input Parameters:
337*8a2c336bSFande Kong +  pc - the HMG context
338*8a2c336bSFande Kong -  reuse - True indicates that HMG will reuse the interpolations
339*8a2c336bSFande Kong 
340*8a2c336bSFande Kong    Options Database Keys:
341*8a2c336bSFande Kong +  -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
342*8a2c336bSFande Kong 
343*8a2c336bSFande Kong    Level: beginner
344*8a2c336bSFande Kong 
345*8a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, reuse, set
346*8a2c336bSFande Kong 
347*8a2c336bSFande Kong .seealso: PCHMG
348*8a2c336bSFande Kong M*/
349*8a2c336bSFande Kong PetscErrorCode PCHMGSetReuseInterpolation(PC pc, PetscBool reuse)
350*8a2c336bSFande Kong {
351*8a2c336bSFande Kong   PC_MG          *mg;
352*8a2c336bSFande Kong   PC_HMG         *hmg;
353*8a2c336bSFande Kong 
354*8a2c336bSFande Kong   PetscFunctionBegin;
355*8a2c336bSFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
356*8a2c336bSFande Kong   mg = (PC_MG*)pc->data;
357*8a2c336bSFande Kong   hmg = (PC_HMG*) mg->innerctx;
358*8a2c336bSFande Kong   hmg->reuseinterp = reuse;
359*8a2c336bSFande Kong   PetscFunctionReturn(0);
360*8a2c336bSFande Kong }
361*8a2c336bSFande Kong 
362*8a2c336bSFande Kong /*MC
363*8a2c336bSFande Kong    PCHMGSetUseSubspaceCoarsening - Use subspace coarsening in HMG
364*8a2c336bSFande Kong 
365*8a2c336bSFande Kong    Logically Collective on PC
366*8a2c336bSFande Kong 
367*8a2c336bSFande Kong    Input Parameters:
368*8a2c336bSFande Kong +  pc - the HMG context
369*8a2c336bSFande Kong -  reuse - True indicates that HMG will use the subspace coarsening
370*8a2c336bSFande Kong 
371*8a2c336bSFande Kong    Options Database Keys:
372*8a2c336bSFande Kong +  -pc_hmg_use_subspace_coarsening  <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix).
373*8a2c336bSFande Kong 
374*8a2c336bSFande Kong    Level: beginner
375*8a2c336bSFande Kong 
376*8a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, subspace, coarsening
377*8a2c336bSFande Kong 
378*8a2c336bSFande Kong .seealso: PCHMG
379*8a2c336bSFande Kong M*/
380*8a2c336bSFande Kong PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC pc, PetscBool subspace)
381*8a2c336bSFande Kong {
382*8a2c336bSFande Kong   PC_MG          *mg;
383*8a2c336bSFande Kong   PC_HMG         *hmg;
384*8a2c336bSFande Kong 
385*8a2c336bSFande Kong   PetscFunctionBegin;
386*8a2c336bSFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
387*8a2c336bSFande Kong   mg = (PC_MG*)pc->data;
388*8a2c336bSFande Kong   hmg = (PC_HMG*) mg->innerctx;
389*8a2c336bSFande Kong   hmg->subcoarsening = subspace;
390*8a2c336bSFande Kong   PetscFunctionReturn(0);
391*8a2c336bSFande Kong }
392*8a2c336bSFande Kong 
393*8a2c336bSFande Kong /*MC
394*8a2c336bSFande Kong    PCHMGSetInnerPCType - Set an inner PC type
395*8a2c336bSFande Kong 
396*8a2c336bSFande Kong    Logically Collective on PC
397*8a2c336bSFande Kong 
398*8a2c336bSFande Kong    Input Parameters:
399*8a2c336bSFande Kong +  pc - the HMG context
400*8a2c336bSFande Kong -  type - <hypre, gamg> coarsening algorithm
401*8a2c336bSFande Kong 
402*8a2c336bSFande Kong    Options Database Keys:
403*8a2c336bSFande Kong +  -pc_hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix
404*8a2c336bSFande Kong 
405*8a2c336bSFande Kong    Level: beginner
406*8a2c336bSFande Kong 
407*8a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, coarsening
408*8a2c336bSFande Kong 
409*8a2c336bSFande Kong .seealso: PCHMG, PCType
410*8a2c336bSFande Kong M*/
411*8a2c336bSFande Kong PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type)
412*8a2c336bSFande Kong {
413*8a2c336bSFande Kong   PC_MG           *mg;
414*8a2c336bSFande Kong   PC_HMG          *hmg;
415*8a2c336bSFande Kong   PetscBool       hypre,gamg;
416*8a2c336bSFande Kong   PetscErrorCode  ierr;
417*8a2c336bSFande Kong 
418*8a2c336bSFande Kong   PetscFunctionBegin;
419*8a2c336bSFande Kong   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
420*8a2c336bSFande Kong   mg = (PC_MG*)pc->data;
421*8a2c336bSFande Kong   hmg = (PC_HMG*) mg->innerctx;
422*8a2c336bSFande Kong   ierr = PetscStrcmp(type,"hypre",&hypre);CHKERRQ(ierr);
423*8a2c336bSFande Kong   ierr = PetscStrcmp(type,"gamg",&gamg);CHKERRQ(ierr);
424*8a2c336bSFande Kong   if (!hypre && !gamg) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Does not support PC type %s \n",type);
425*8a2c336bSFande Kong   ierr = PetscStrallocpy(type,&(hmg->innerpctype));CHKERRQ(ierr);
426*8a2c336bSFande Kong   PetscFunctionReturn(0);
427*8a2c336bSFande Kong }
428*8a2c336bSFande Kong 
429*8a2c336bSFande Kong /*MC
430360ee056SFande Kong    PCHMG - Hybrid of PETSc preconditioners (such as ASM, BJacobi, SOR, etc.) and Hypre BoomerAMG. BoomerAMG is used to coarsen matrix and generate
431360ee056SFande Kong            a sequence of coarse matrices and interpolations. The matrices and interpolations are employed to construct PCMG, and then any available
432360ee056SFande Kong            PETSc preconditioners can be chosen as smoothers and the coarse solver.
433360ee056SFande Kong 
434360ee056SFande Kong    Options Database Keys:
435*8a2c336bSFande Kong +  -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
436*8a2c336bSFande Kong .  -pc_hmg_use_subspace_coarsening  <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix).
437*8a2c336bSFande Kong .  -pc_hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix
438360ee056SFande Kong 
439360ee056SFande Kong 
440360ee056SFande Kong    Notes:
441360ee056SFande Kong     For multicomponent problems, we can just coarsen one submatrix associated with one particular component. In this way, the preconditioner setup
442360ee056SFande Kong     time is significantly reduced. One typical use case is neutron transport equations. There are many variables on each mesh vertex due to the
443360ee056SFande Kong     of angle and energy. Each variable, in fact, corresponds to the same PDEs but with different material properties.
444360ee056SFande Kong 
445360ee056SFande Kong    Level: beginner
446360ee056SFande Kong 
447*8a2c336bSFande Kong    Concepts: Hybrid of ASM and MG, Subspace Coarsening
448360ee056SFande Kong 
449360ee056SFande Kong     References:
450360ee056SFande 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
451360ee056SFande Kong     Newton-Krylov-Schwarz method with subspace-based coarsening and partition-based balancing for the multigroup neutron transport equations on
452360ee056SFande Kong     3D unstructured meshes, arXiv preprint arXiv:1903.03659, 2019
453360ee056SFande Kong 
454*8a2c336bSFande Kong .seealso:  PCCreate(), PCSetType(), PCType, PC, PCMG, PCHYPRE, PCHMG, PCHYPREBoomerAMGGetCoarseOperators(), PCHYPREBoomerAMGGetInterpolations(), PCMGGetInterpolations(), PCMGGetCoarseOperators()
455360ee056SFande Kong 
456360ee056SFande Kong M*/
457360ee056SFande Kong PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc)
458360ee056SFande Kong {
459360ee056SFande Kong   PetscErrorCode ierr;
460360ee056SFande Kong   PC_HMG         *hmg;
461360ee056SFande Kong   PC_MG          *mg;
462360ee056SFande Kong 
463360ee056SFande Kong   PetscFunctionBegin;
464360ee056SFande Kong   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
465360ee056SFande Kong   if (pc->ops->destroy) {
466360ee056SFande Kong     ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr);
467360ee056SFande Kong     pc->data = 0;
468360ee056SFande Kong   }
469360ee056SFande Kong   ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr);
470360ee056SFande Kong 
471360ee056SFande Kong   ierr         = PCSetType(pc,PCMG);CHKERRQ(ierr);
472360ee056SFande Kong   ierr         = PetscNew(&hmg);CHKERRQ(ierr); \
473360ee056SFande Kong 
474360ee056SFande Kong   mg                      = (PC_MG*) pc->data;
475360ee056SFande Kong   mg->innerctx            = hmg;
476360ee056SFande Kong   hmg->reuseinterp        = PETSC_FALSE;
477*8a2c336bSFande Kong   hmg->subcoarsening      = PETSC_FALSE;
478*8a2c336bSFande Kong   hmg->innerpc            = NULL;
479360ee056SFande Kong 
480360ee056SFande Kong   pc->ops->setfromoptions = PCSetFromOptions_HMG;
481360ee056SFande Kong   pc->ops->view           = PCView_HMG;
482360ee056SFande Kong   pc->ops->destroy        = PCDestroy_HMG;
483360ee056SFande Kong   pc->ops->setup          = PCSetUp_HMG;
484360ee056SFande Kong   PetscFunctionReturn(0);
485360ee056SFande Kong }
486