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