xref: /petsc/src/ksp/pc/impls/hmg/hmg.c (revision b74e562f2c2fccefb9b214c8577f9611b088b129)
1 #include <petscdm.h>
2 #include <petscctable.h>
3 #include <petsc/private/matimpl.h>
4 #include <petsc/private/pcmgimpl.h>
5 #include <petsc/private/pcimpl.h>      /*I "petscpc.h" I*/
6 
7 typedef struct {
8   PC               innerpc;            /* A MG inner PC (Hypre or PCGAMG) to setup interpolations and coarse operators  */
9   char*            innerpctype;        /* PCGAMG or PCHYPRE */
10   PetscBool        reuseinterp;        /* A flag indicates if or not to reuse the interpolations */
11   PetscBool        subcoarsening;
12   PetscErrorCode   (*getoperators)(PC,PetscInt*,Mat**);
13   PetscErrorCode   (*getinterpolations)(PC,PetscInt*,Mat**);
14   PetscErrorCode   (*setfromoptions)(PetscOptionItems*,PC);
15 } PC_HMG;
16 
17 PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems*,PC);
18 PetscErrorCode PCReset_MG(PC);
19 #if PETSC_HAVE_HYPRE
20 PetscErrorCode PCHYPREBoomerAMGGetCoarseOperators(PC,PetscInt*,Mat**);
21 PetscErrorCode PCHYPREBoomerAMGGetInterpolations(PC,PetscInt*,Mat**);
22 PetscErrorCode PCSetFromOptions_HYPRE(PetscOptionItems*,PC);
23 #endif
24 PetscErrorCode PCMGGetCoarseOperators(PC,PetscInt*,Mat**);
25 PetscErrorCode PCMGGetInterpolations(PC,PetscInt*,Mat **);
26 PetscErrorCode PCSetFromOptions_GAMG(PetscOptionItems*,PC);
27 
28 static PetscErrorCode PCHMGSetInnerPCFromOptions_Private(PC pc,PetscOptionItems *PetscOptionsObject)
29 {
30   PetscErrorCode     ierr;
31   PC_MG              *mg   = (PC_MG*)pc->data;
32   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
33 
34   PetscFunctionBegin;
35   if (hmg->setfromoptions) {
36     ierr = (*hmg->setfromoptions)(PetscOptionsObject,hmg->innerpc);CHKERRQ(ierr);
37   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_NULL,"Inner PC does not support SetFromOptions\n");
38   PetscFunctionReturn(0);
39 }
40 
41 static PetscErrorCode PCHMGGetInnerPCCoarseOperators_Private(PC pc,PetscInt *nlevels,Mat **operators)
42 {
43   PetscErrorCode     ierr;
44   PC_MG              *mg   = (PC_MG*)pc->data;
45   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
46 
47   PetscFunctionBegin;
48   if (hmg->getoperators) {
49     ierr = (*hmg->getoperators)(hmg->innerpc,nlevels,operators);CHKERRQ(ierr);
50   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_NULL,"Inner PC does not support GetOperators\n");
51   PetscFunctionReturn(0);
52 }
53 
54 static PetscErrorCode PCHMGGetInnerPCInterpolations_Private(PC pc,PetscInt *nlevels,Mat *interpolations[])
55 {
56   PetscErrorCode     ierr;
57   PC_MG              *mg   = (PC_MG*)pc->data;
58   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
59 
60   PetscFunctionBegin;
61   if (hmg->getinterpolations) {
62     ierr = (*hmg->getinterpolations)(hmg->innerpc,nlevels,interpolations);CHKERRQ(ierr);
63   } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_NULL,"Inner PC does not support GetInterpolation\n");
64   PetscFunctionReturn(0);
65 }
66 
67 static PetscErrorCode PCHMGExtractSubMatrix_Private(Mat pmat,Mat *submat,MatReuse reuse,PetscInt blocksize)
68 {
69   IS             isrow;
70   PetscErrorCode ierr;
71   PetscInt       rstart,rend;
72   MPI_Comm       comm;
73 
74   PetscFunctionBegin;
75   ierr = PetscObjectGetComm((PetscObject)pmat,&comm);CHKERRQ(ierr);
76   ierr = MatGetOwnershipRange(pmat,&rstart,&rend);CHKERRQ(ierr);
77   if ((rend-rstart)%blocksize != 0) SETERRQ3(comm,PETSC_ERR_ARG_INCOMP,"Block size %d is inconsisent for [%d, %d) \n",blocksize,rstart,rend);
78   ierr = ISCreateStride(comm,(rend-rstart)/blocksize,rstart,blocksize,&isrow);
79   ierr = MatCreateSubMatrix(pmat,isrow,isrow,reuse,submat);CHKERRQ(ierr);
80   ierr = ISDestroy(&isrow);CHKERRQ(ierr);
81   PetscFunctionReturn(0);
82 }
83 
84 static PetscErrorCode PCHMGExpandInterpolation_Private(Mat subinterp, Mat *interp, PetscInt blocksize)
85 {
86   PetscInt              subrstart,subrend,subrowsize,subcolsize,subcstart,subcend,rowsize,colsize;
87   PetscInt              subrow,row,nz,*d_nnz,*o_nnz,i,j,dnz,onz,max_nz,*indices;
88   const PetscInt        *idx;
89   const PetscScalar     *values;
90   PetscErrorCode        ierr;
91   MPI_Comm              comm;
92 
93   PetscFunctionBegin;
94   ierr = PetscObjectGetComm((PetscObject)subinterp,&comm);CHKERRQ(ierr);
95   ierr = MatGetOwnershipRange(subinterp,&subrstart,&subrend);CHKERRQ(ierr);
96   subrowsize = subrend-subrstart;
97   rowsize = subrowsize*blocksize;
98   ierr = PetscCalloc2(rowsize,&d_nnz,rowsize,&o_nnz);CHKERRQ(ierr);
99   ierr = MatGetOwnershipRangeColumn(subinterp,&subcstart,&subcend);CHKERRQ(ierr);
100   subcolsize = subcend - subcstart;
101   colsize    = subcolsize*blocksize;
102   max_nz = 0;
103   for (subrow=subrstart;subrow<subrend;subrow++) {
104     ierr = MatGetRow(subinterp,subrow,&nz,&idx,NULL);CHKERRQ(ierr);
105     if (max_nz<nz) max_nz = nz;
106     dnz = 0; onz = 0;
107     for (i=0;i<nz;i++) {
108       if(idx[i]>=subcstart && idx[i]<subcend) dnz++;
109       else onz++;
110     }
111     for (i=0;i<blocksize;i++) {
112       d_nnz[(subrow-subrstart)*blocksize+i] = dnz;
113       o_nnz[(subrow-subrstart)*blocksize+i] = onz;
114     }
115     ierr = MatRestoreRow(subinterp,subrow,&nz,&idx,NULL);CHKERRQ(ierr);
116   }
117   ierr = MatCreateAIJ(comm,rowsize,colsize,PETSC_DETERMINE,PETSC_DETERMINE,0,d_nnz,0,o_nnz,interp);CHKERRQ(ierr);
118   ierr = MatSetOption(*interp,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
119   ierr = MatSetOption(*interp,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr);
120   ierr = MatSetOption(*interp,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
121   ierr = MatSetFromOptions(*interp);CHKERRQ(ierr);
122 
123   ierr = MatSetUp(*interp);CHKERRQ(ierr);
124   ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr);
125   ierr = PetscCalloc1(max_nz,&indices);CHKERRQ(ierr);
126   for (subrow=subrstart; subrow<subrend; subrow++) {
127     ierr = MatGetRow(subinterp,subrow,&nz,&idx,&values);CHKERRQ(ierr);
128     for (i=0;i<blocksize;i++) {
129       row = subrow*blocksize+i;
130       for (j=0;j<nz;j++) {
131         indices[j] = idx[j]*blocksize+i;
132       }
133       ierr = MatSetValues(*interp,1,&row,nz,indices,values,INSERT_VALUES);CHKERRQ(ierr);
134     }
135     ierr = MatRestoreRow(subinterp,subrow,&nz,&idx,&values);CHKERRQ(ierr);
136   }
137   ierr = PetscFree(indices);CHKERRQ(ierr);
138   ierr = MatAssemblyBegin(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139   ierr = MatAssemblyEnd(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
140   PetscFunctionReturn(0);
141 }
142 
143 PetscErrorCode PCSetUp_HMG(PC pc)
144 {
145   PetscErrorCode     ierr;
146   Mat                PA, submat;
147   PC_MG              *mg   = (PC_MG*)pc->data;
148   PC_HMG             *hmg   = (PC_HMG*) mg->innerctx;
149   MPI_Comm           comm;
150   PetscInt           level;
151   PetscInt           num_levels;
152   PetscBool          same;
153   Mat                *operators,*interpolations;
154   PetscInt           blocksize;
155 
156   PetscFunctionBegin;
157   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
158   if (pc->setupcalled) {
159    /* Only reuse interpolations when nonzero pattern does not change */
160    if (pc->flag == SAME_NONZERO_PATTERN && hmg->reuseinterp) {
161      ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT);CHKERRQ(ierr);
162      ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
163      PetscFunctionReturn(0);
164     } else {
165      ierr = PCReset_MG(pc);CHKERRQ(ierr);
166      pc->setupcalled = PETSC_FALSE;
167     }
168   }
169 
170   /* Create an inner PC (GAMG or HYPRE) */
171   if (!hmg->innerpc) {
172     ierr = PCCreate(comm,&hmg->innerpc);CHKERRQ(ierr);
173     ierr = PCSetType(hmg->innerpc,hmg->innerpctype);CHKERRQ(ierr);
174     ierr = PetscStrcmp(hmg->innerpctype,"hypre",&same);CHKERRQ(ierr);
175     if (same) {
176       ierr = PCHYPRESetType(hmg->innerpc,"boomeramg");CHKERRQ(ierr);
177     }
178   }
179   ierr = PCGetOperators(pc,NULL,&PA);CHKERRQ(ierr);
180   /* Users need to correctly set a block size of matrix in order to use subspace coarsening */
181   ierr = MatGetBlockSize(PA,&blocksize);CHKERRQ(ierr);
182   if (blocksize<=1) hmg->subcoarsening = PETSC_FALSE;
183   /* Extract a submatrix for constructing subinterpolations */
184   if (hmg->subcoarsening) {
185     ierr = PCHMGExtractSubMatrix_Private(PA,&submat,MAT_INITIAL_MATRIX,blocksize);CHKERRQ(ierr);
186     PA = submat;
187   }
188   ierr = PCSetOperators(hmg->innerpc,PA,PA);CHKERRQ(ierr);
189   if (hmg->subcoarsening) {
190    ierr = MatDestroy(&PA);CHKERRQ(ierr);
191   }
192   /* Setup inner PC correctly. During this step, matrix will be coarsened */
193   ierr = PCSetUseAmat(hmg->innerpc,PETSC_FALSE);CHKERRQ(ierr);
194   ierr = PetscObjectOptionsBegin((PetscObject)hmg->innerpc);CHKERRQ(ierr);
195   ierr = PCHMGSetInnerPCFromOptions_Private(pc,PetscOptionsObject);CHKERRQ(ierr);
196   ierr = PetscOptionsEnd();CHKERRQ(ierr);
197   ierr = PCSetUp(hmg->innerpc);
198 
199   /* Obtain interpolations IN PLACE. For BoomerAMG, (I,J,data) is reused to avoid memory overhead */
200   ierr = PCHMGGetInnerPCInterpolations_Private(pc,&num_levels,&interpolations);CHKERRQ(ierr);
201   /* We can reuse the coarse operators when we do the full space coarsening */
202   if (!hmg->subcoarsening) {
203     ierr = PCHMGGetInnerPCCoarseOperators_Private(pc,&num_levels,&operators);CHKERRQ(ierr);
204   }
205 
206   ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr);
207   hmg->innerpc = 0;
208   ierr = PCMGSetLevels_MG(pc,num_levels,NULL);CHKERRQ(ierr);
209   /* Set coarse matrices and interpolations to PCMG */
210   for (level=num_levels-1; level>0; level--) {
211     Mat P=0, pmat=0;
212     Vec b, x,r;
213     if (hmg->subcoarsening) {
214       /* Grow interpolation. In the future, we should use MAIJ */
215       ierr = PCHMGExpandInterpolation_Private(interpolations[level-1],&P,blocksize);CHKERRQ(ierr);
216       ierr = MatDestroy(&interpolations[level-1]);CHKERRQ(ierr);
217     } else {
218       P = interpolations[level-1];
219     }
220     ierr = MatCreateVecs(P,&b,&r);CHKERRQ(ierr);
221     ierr = PCMGSetInterpolation(pc,level,P);CHKERRQ(ierr);
222     ierr = PCMGSetRestriction(pc,level,P);CHKERRQ(ierr);
223     ierr = MatDestroy(&P);CHKERRQ(ierr);
224     /* We reuse the matrices when we do not do subspace coarsening */
225     if ((level-1)>=0 && !hmg->subcoarsening) {
226       pmat = operators[level-1];
227       ierr = PCMGSetOperators(pc,level-1,pmat,pmat);CHKERRQ(ierr);
228       ierr = MatDestroy(&pmat);CHKERRQ(ierr);
229     }
230     ierr = PCMGSetRhs(pc,level-1,b);CHKERRQ(ierr);
231 
232     ierr = PCMGSetR(pc,level,r);CHKERRQ(ierr);
233     ierr = VecDestroy(&r);CHKERRQ(ierr);
234 
235     ierr = VecDuplicate(b,&x);
236     ierr = PCMGSetX(pc,level-1,x);CHKERRQ(ierr);
237     ierr = VecDestroy(&x);CHKERRQ(ierr);
238     ierr = VecDestroy(&b);CHKERRQ(ierr);
239   }
240   ierr = PetscFree(interpolations);CHKERRQ(ierr);
241   if (!hmg->subcoarsening) {
242     ierr = PetscFree(operators);CHKERRQ(ierr);
243   }
244   /* Turn Galerkin off when we already have coarse operators */
245   ierr = PCMGSetGalerkin(pc,hmg->subcoarsening ? PC_MG_GALERKIN_PMAT:PC_MG_GALERKIN_NONE);CHKERRQ(ierr);
246   ierr = PCSetDM(pc,NULL);CHKERRQ(ierr);
247   ierr = PCSetUseAmat(pc,PETSC_FALSE);CHKERRQ(ierr);
248   ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr);
249   ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */
250   ierr = PetscOptionsEnd();CHKERRQ(ierr);
251   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
252   PetscFunctionReturn(0);
253 }
254 
255 PetscErrorCode PCDestroy_HMG(PC pc)
256 {
257   PetscErrorCode ierr;
258   PC_MG          *mg  = (PC_MG*)pc->data;
259   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
260 
261   PetscFunctionBegin;
262   ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr);
263   ierr = PetscFree(hmg->innerpctype);CHKERRQ(ierr);
264   ierr = PetscFree(hmg);CHKERRQ(ierr);
265   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
266   PetscFunctionReturn(0);
267 }
268 
269 PetscErrorCode PCView_HMG(PC pc,PetscViewer viewer)
270 {
271   PC_MG          *mg = (PC_MG*)pc->data;
272   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
273   PetscErrorCode ierr;
274   PetscBool      iascii;
275 
276   PetscFunctionBegin;
277   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
278   if (iascii) {
279     ierr = PetscViewerASCIIPrintf(viewer," Reuse interpolation: %s\n",hmg->reuseinterp? "true":"false");CHKERRQ(ierr);
280     ierr = PetscViewerASCIIPrintf(viewer," Use subspace coarsening: %s\n",hmg->subcoarsening? "true":"false");CHKERRQ(ierr);
281     ierr = PetscViewerASCIIPrintf(viewer," Inner PC type: %s \n",hmg->innerpctype);CHKERRQ(ierr);
282   }
283   ierr = PCView_MG(pc,viewer);CHKERRQ(ierr);
284   PetscFunctionReturn(0);
285 }
286 
287 #if PETSC_HAVE_HYPRE
288 static const char *InnerPCType[]   = {"gamg","hypre"};
289 #else
290 static const char *InnerPCType[]   = {"gamg"};
291 #endif
292 
293 PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems *PetscOptionsObject,PC pc)
294 {
295   PetscErrorCode ierr;
296   PC_MG          *mg = (PC_MG*)pc->data;
297   PC_HMG         *hmg = (PC_HMG*) mg->innerctx;
298   PetscBool      flg;
299   PetscInt       indx;
300 
301   PetscFunctionBegin;
302   ierr = PetscOptionsHead(PetscOptionsObject,"HMG");CHKERRQ(ierr);
303   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   ierr = PetscOptionsBool("-pc_hmg_use_subspace_coarsening","Use the subspace coarsening to compute the interpolations","None",hmg->subcoarsening,&hmg->subcoarsening,NULL);CHKERRQ(ierr);
305 #if PETSC_HAVE_HYPRE
306   indx = 1;
307 #else
308   indx = 0;
309 #endif
310   ierr = PetscOptionsEList("-pc_hmg_inner_pc_type","Inner PC type",NULL,InnerPCType,2,"gamg",&indx,&flg);CHKERRQ(ierr);
311  if (indx == 0) {
312    hmg->getinterpolations = PCMGGetInterpolations;
313    hmg->getoperators      = PCMGGetCoarseOperators;
314    hmg->setfromoptions    = PCSetFromOptions_GAMG;
315    ierr = PetscStrallocpy(PCGAMG,&(hmg->innerpctype));CHKERRQ(ierr);
316   }
317 #if PETSC_HAVE_HYPRE
318  else if (indx == 1) {
319     hmg->getinterpolations = PCHYPREBoomerAMGGetInterpolations;
320     hmg->getoperators      = PCHYPREBoomerAMGGetCoarseOperators;
321     hmg->setfromoptions    = PCSetFromOptions_HYPRE;
322     ierr = PetscStrallocpy(PCHYPRE,&(hmg->innerpctype));CHKERRQ(ierr);
323   }
324 #endif
325   else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown inner PC type");
326 
327   ierr = PetscOptionsTail();CHKERRQ(ierr);
328   PetscFunctionReturn(0);
329 }
330 
331 /*MC
332    PCHMGSetReuseInterpolation - Reuse interpolation matrices in HMG
333 
334    Logically Collective on PC
335 
336    Input Parameters:
337 +  pc - the HMG context
338 -  reuse - True indicates that HMG will reuse the interpolations
339 
340    Options Database Keys:
341 +  -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
342 
343    Level: beginner
344 
345 .keywords: HMG, multigrid, interpolation, reuse, set
346 
347 .seealso: PCHMG
348 M*/
349 PetscErrorCode PCHMGSetReuseInterpolation(PC pc, PetscBool reuse)
350 {
351   PC_MG          *mg;
352   PC_HMG         *hmg;
353 
354   PetscFunctionBegin;
355   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
356   mg = (PC_MG*)pc->data;
357   hmg = (PC_HMG*) mg->innerctx;
358   hmg->reuseinterp = reuse;
359   PetscFunctionReturn(0);
360 }
361 
362 /*MC
363    PCHMGSetUseSubspaceCoarsening - Use subspace coarsening in HMG
364 
365    Logically Collective on PC
366 
367    Input Parameters:
368 +  pc - the HMG context
369 -  reuse - True indicates that HMG will use the subspace coarsening
370 
371    Options Database Keys:
372 +  -pc_hmg_use_subspace_coarsening  <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix).
373 
374    Level: beginner
375 
376 .keywords: HMG, multigrid, interpolation, subspace, coarsening
377 
378 .seealso: PCHMG
379 M*/
380 PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC pc, PetscBool subspace)
381 {
382   PC_MG          *mg;
383   PC_HMG         *hmg;
384 
385   PetscFunctionBegin;
386   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
387   mg = (PC_MG*)pc->data;
388   hmg = (PC_HMG*) mg->innerctx;
389   hmg->subcoarsening = subspace;
390   PetscFunctionReturn(0);
391 }
392 
393 /*MC
394    PCHMGSetInnerPCType - Set an inner PC type
395 
396    Logically Collective on PC
397 
398    Input Parameters:
399 +  pc - the HMG context
400 -  type - <hypre, gamg> coarsening algorithm
401 
402    Options Database Keys:
403 +  -pc_hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix
404 
405    Level: beginner
406 
407 .keywords: HMG, multigrid, interpolation, coarsening
408 
409 .seealso: PCHMG, PCType
410 M*/
411 PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type)
412 {
413   PC_MG           *mg;
414   PC_HMG          *hmg;
415   PetscBool       hypre,gamg;
416   PetscErrorCode  ierr;
417 
418   PetscFunctionBegin;
419   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
420   mg = (PC_MG*)pc->data;
421   hmg = (PC_HMG*) mg->innerctx;
422   ierr = PetscStrcmp(type,"hypre",&hypre);CHKERRQ(ierr);
423   ierr = PetscStrcmp(type,"gamg",&gamg);CHKERRQ(ierr);
424   if (!hypre && !gamg) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_UNKNOWN_TYPE,"Does not support PC type %s \n",type);
425   ierr = PetscStrallocpy(type,&(hmg->innerpctype));CHKERRQ(ierr);
426   PetscFunctionReturn(0);
427 }
428 
429 /*MC
430    PCHMG - Hybrid of PETSc preconditioners (such as ASM, BJacobi, SOR, etc.) and Hypre BoomerAMG. BoomerAMG is used to coarsen matrix and generate
431            a sequence of coarse matrices and interpolations. The matrices and interpolations are employed to construct PCMG, and then any available
432            PETSc preconditioners can be chosen as smoothers and the coarse solver.
433 
434    Options Database Keys:
435 +  -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time.
436 .  -pc_hmg_use_subspace_coarsening  <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix).
437 .  -pc_hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix
438 
439 
440    Notes:
441     For multicomponent problems, we can just coarsen one submatrix associated with one particular component. In this way, the preconditioner setup
442     time is significantly reduced. One typical use case is neutron transport equations. There are many variables on each mesh vertex due to the
443     of angle and energy. Each variable, in fact, corresponds to the same PDEs but with different material properties.
444 
445    Level: beginner
446 
447    Concepts: Hybrid of ASM and MG, Subspace Coarsening
448 
449     References:
450 +   1. - Fande Kong, Yaqi Wang, Derek R Gaston, Cody J Permann, Andrew E Slaughter, Alexander D Lindsay, Richard C Martineau, A highly parallel multilevel
451     Newton-Krylov-Schwarz method with subspace-based coarsening and partition-based balancing for the multigroup neutron transport equations on
452     3D unstructured meshes, arXiv preprint arXiv:1903.03659, 2019
453 
454 .seealso:  PCCreate(), PCSetType(), PCType, PC, PCMG, PCHYPRE, PCHMG, PCHYPREBoomerAMGGetCoarseOperators(), PCHYPREBoomerAMGGetInterpolations(), PCMGGetInterpolations(), PCMGGetCoarseOperators()
455 
456 M*/
457 PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc)
458 {
459   PetscErrorCode ierr;
460   PC_HMG         *hmg;
461   PC_MG          *mg;
462 
463   PetscFunctionBegin;
464   /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */
465   if (pc->ops->destroy) {
466     ierr =  (*pc->ops->destroy)(pc);CHKERRQ(ierr);
467     pc->data = 0;
468   }
469   ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr);
470 
471   ierr         = PCSetType(pc,PCMG);CHKERRQ(ierr);
472   ierr         = PetscNew(&hmg);CHKERRQ(ierr); \
473 
474   mg                      = (PC_MG*) pc->data;
475   mg->innerctx            = hmg;
476   hmg->reuseinterp        = PETSC_FALSE;
477   hmg->subcoarsening      = PETSC_FALSE;
478   hmg->innerpc            = NULL;
479 
480   pc->ops->setfromoptions = PCSetFromOptions_HMG;
481   pc->ops->view           = PCView_HMG;
482   pc->ops->destroy        = PCDestroy_HMG;
483   pc->ops->setup          = PCSetUp_HMG;
484   PetscFunctionReturn(0);
485 }
486