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 PCMGGalerkinType galerkin; 106 107 PetscFunctionBegin; 108 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 109 if (pc->setupcalled) { 110 /* Only reuse interpolations when nonzero pattern does not change 111 * since coarse matrices are already allocated. 112 * */ 113 if (pc->flag == SAME_NONZERO_PATTERN && hmg->reuseinterp) { 114 /* If we did not use Galerkin in the last call, we have to build from scratch */ 115 ierr = PCMGGetGalerkin(pc,&galerkin);CHKERRQ(ierr); 116 if (galerkin == PC_MG_GALERKIN_NONE) 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 /* Grow interpolation. In the future, we should use MAIJ */ 177 ierr = PCHMGExpandInterpolation_Private(interpolations[level-1],&P,blocksize);CHKERRQ(ierr); 178 ierr = MatDestroy(&interpolations[level-1]);CHKERRQ(ierr); 179 } else { 180 P = interpolations[level-1]; 181 } 182 ierr = MatCreateVecs(P,&b,&r);CHKERRQ(ierr); 183 ierr = PCMGSetInterpolation(pc,level,P);CHKERRQ(ierr); 184 ierr = PCMGSetRestriction(pc,level,P);CHKERRQ(ierr); 185 ierr = MatDestroy(&P);CHKERRQ(ierr); 186 /* We reuse the matrices when we do not do subspace coarsening */ 187 if ((level-1)>=0 && !hmg->subcoarsening) { 188 pmat = operators[level-1]; 189 ierr = PCMGSetOperators(pc,level-1,pmat,pmat);CHKERRQ(ierr); 190 ierr = MatDestroy(&pmat);CHKERRQ(ierr); 191 } 192 ierr = PCMGSetRhs(pc,level-1,b);CHKERRQ(ierr); 193 194 ierr = PCMGSetR(pc,level,r);CHKERRQ(ierr); 195 ierr = VecDestroy(&r);CHKERRQ(ierr); 196 197 ierr = VecDuplicate(b,&x);CHKERRQ(ierr); 198 ierr = PCMGSetX(pc,level-1,x);CHKERRQ(ierr); 199 ierr = VecDestroy(&x);CHKERRQ(ierr); 200 ierr = VecDestroy(&b);CHKERRQ(ierr); 201 } 202 ierr = PetscFree(interpolations);CHKERRQ(ierr); 203 if (!hmg->subcoarsening) { 204 ierr = PetscFree(operators);CHKERRQ(ierr); 205 } 206 /* Turn Galerkin off when we already have coarse operators */ 207 ierr = PCMGSetGalerkin(pc,hmg->subcoarsening ? PC_MG_GALERKIN_PMAT:PC_MG_GALERKIN_NONE);CHKERRQ(ierr); 208 ierr = PCSetDM(pc,NULL);CHKERRQ(ierr); 209 ierr = PCSetUseAmat(pc,PETSC_FALSE);CHKERRQ(ierr); 210 ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 211 ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */ 212 ierr = PetscOptionsEnd();CHKERRQ(ierr); 213 ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 214 PetscFunctionReturn(0); 215 } 216 217 PetscErrorCode PCDestroy_HMG(PC pc) 218 { 219 PetscErrorCode ierr; 220 PC_MG *mg = (PC_MG*)pc->data; 221 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 222 223 PetscFunctionBegin; 224 ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr); 225 ierr = PetscFree(hmg->innerpctype);CHKERRQ(ierr); 226 ierr = PetscFree(hmg);CHKERRQ(ierr); 227 ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 228 229 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",NULL);CHKERRQ(ierr); 230 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",NULL);CHKERRQ(ierr); 231 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",NULL);CHKERRQ(ierr); 232 PetscFunctionReturn(0); 233 } 234 235 PetscErrorCode PCView_HMG(PC pc,PetscViewer viewer) 236 { 237 PC_MG *mg = (PC_MG*)pc->data; 238 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 239 PetscErrorCode ierr; 240 PetscBool iascii; 241 242 PetscFunctionBegin; 243 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 244 if (iascii) { 245 ierr = PetscViewerASCIIPrintf(viewer," Reuse interpolation: %s\n",hmg->reuseinterp? "true":"false");CHKERRQ(ierr); 246 ierr = PetscViewerASCIIPrintf(viewer," Use subspace coarsening: %s\n",hmg->subcoarsening? "true":"false");CHKERRQ(ierr); 247 ierr = PetscViewerASCIIPrintf(viewer," Inner PC type: %s \n",hmg->innerpctype);CHKERRQ(ierr); 248 } 249 ierr = PCView_MG(pc,viewer);CHKERRQ(ierr); 250 PetscFunctionReturn(0); 251 } 252 253 PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems *PetscOptionsObject,PC pc) 254 { 255 PetscErrorCode ierr; 256 PC_MG *mg = (PC_MG*)pc->data; 257 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 258 259 PetscFunctionBegin; 260 ierr = PetscOptionsHead(PetscOptionsObject,"HMG");CHKERRQ(ierr); 261 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); 262 ierr = PetscOptionsBool("-pc_hmg_use_subspace_coarsening","Use the subspace coarsening to compute the interpolations","None",hmg->subcoarsening,&hmg->subcoarsening,NULL);CHKERRQ(ierr); 263 ierr = PetscOptionsTail();CHKERRQ(ierr); 264 PetscFunctionReturn(0); 265 } 266 267 static PetscErrorCode PCHMGSetReuseInterpolation_HMG(PC pc, PetscBool reuse) 268 { 269 PC_MG *mg = (PC_MG*)pc->data; 270 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 271 272 PetscFunctionBegin; 273 hmg->reuseinterp = reuse; 274 PetscFunctionReturn(0); 275 } 276 277 /*MC 278 PCHMGSetReuseInterpolation - Reuse interpolation matrices in HMG 279 280 Logically Collective on PC 281 282 Input Parameters: 283 + pc - the HMG context 284 - reuse - True indicates that HMG will reuse the interpolations 285 286 Options Database Keys: 287 + -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time. 288 289 Level: beginner 290 291 .keywords: HMG, multigrid, interpolation, reuse, set 292 293 .seealso: PCHMG 294 M*/ 295 PetscErrorCode PCHMGSetReuseInterpolation(PC pc, PetscBool reuse) 296 { 297 PetscErrorCode ierr; 298 299 PetscFunctionBegin; 300 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 301 ierr = PetscUseMethod(pc,"PCHMGSetReuseInterpolation_C",(PC,PetscBool),(pc,reuse));CHKERRQ(ierr); 302 PetscFunctionReturn(0); 303 } 304 305 static PetscErrorCode PCHMGSetUseSubspaceCoarsening_HMG(PC pc, PetscBool subspace) 306 { 307 PC_MG *mg = (PC_MG*)pc->data; 308 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 309 310 PetscFunctionBegin; 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 = (PC_MG*)pc->data; 346 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 347 PetscErrorCode ierr; 348 349 PetscFunctionBegin; 350 ierr = PetscStrallocpy(type,&(hmg->innerpctype));CHKERRQ(ierr); 351 PetscFunctionReturn(0); 352 } 353 354 /*MC 355 PCHMGSetInnerPCType - Set an inner PC type 356 357 Logically Collective on PC 358 359 Input Parameters: 360 + pc - the HMG context 361 - type - <hypre, gamg> coarsening algorithm 362 363 Options Database Keys: 364 + -hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix 365 366 Level: beginner 367 368 .keywords: HMG, multigrid, interpolation, coarsening 369 370 .seealso: PCHMG, PCType 371 M*/ 372 PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type) 373 { 374 PetscErrorCode ierr; 375 376 PetscFunctionBegin; 377 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 378 ierr = PetscUseMethod(pc,"PCHMGSetInnerPCType_C",(PC,PCType),(pc,type));CHKERRQ(ierr); 379 PetscFunctionReturn(0); 380 } 381 382 /*MC 383 PCHMG - Hybrid of PETSc preconditioners (such as ASM, BJacobi, SOR, etc.) and Hypre BoomerAMG, GAMG or other multilevel methods. BoomerAMG, GAMG 384 or other multilevel methods is used to coarsen matrix and generate a sequence of coarse matrices and interpolations. The matrices and 385 interpolations are employed to construct PCMG, and then any available PETSc preconditioners can be chosen as smoothers and the coarse solver. 386 387 Options Database Keys: 388 + -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time. 389 . -pc_hmg_use_subspace_coarsening <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix). 390 . -hmg_inner_pc_type <hypre, gamg, ...> - What method is used to coarsen matrix 391 392 393 Notes: 394 For multicomponent problems, we can just coarsen one submatrix associated with one particular component. In this way, the preconditioner setup 395 time is significantly reduced. One typical use case is neutron transport equations. There are many variables on each mesh vertex due to the 396 of angle and energy. Each variable, in fact, corresponds to the same PDEs but with different material properties. 397 398 Level: beginner 399 400 Concepts: Hybrid of ASM and MG, Subspace Coarsening 401 402 References: 403 + 1. - Fande Kong, Yaqi Wang, Derek R Gaston, Cody J Permann, Andrew E Slaughter, Alexander D Lindsay, Richard C Martineau, A highly parallel multilevel 404 Newton-Krylov-Schwarz method with subspace-based coarsening and partition-based balancing for the multigroup neutron transport equations on 405 3D unstructured meshes, arXiv preprint arXiv:1903.03659, 2019 406 407 .seealso: PCCreate(), PCSetType(), PCType, PC, PCMG, PCHYPRE, PCHMG, PCGetCoarseOperators(), PCGetInterpolations(), PCHMGSetReuseInterpolation(), PCHMGSetUseSubspaceCoarsening(), 408 PCHMGSetInnerPCType() 409 410 M*/ 411 PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc) 412 { 413 PetscErrorCode ierr; 414 PC_HMG *hmg; 415 PC_MG *mg; 416 417 PetscFunctionBegin; 418 /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */ 419 if (pc->ops->destroy) { 420 ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 421 pc->data = 0; 422 } 423 ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr); 424 425 ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); 426 ierr = PetscNew(&hmg);CHKERRQ(ierr); 427 428 mg = (PC_MG*) pc->data; 429 mg->innerctx = hmg; 430 hmg->reuseinterp = PETSC_FALSE; 431 hmg->subcoarsening = PETSC_FALSE; 432 hmg->innerpc = NULL; 433 434 pc->ops->setfromoptions = PCSetFromOptions_HMG; 435 pc->ops->view = PCView_HMG; 436 pc->ops->destroy = PCDestroy_HMG; 437 pc->ops->setup = PCSetUp_HMG; 438 439 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",PCHMGSetReuseInterpolation_HMG);CHKERRQ(ierr); 440 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",PCHMGSetUseSubspaceCoarsening_HMG);CHKERRQ(ierr); 441 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",PCHMGSetInnerPCType_HMG);CHKERRQ(ierr); 442 PetscFunctionReturn(0); 443 } 444