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