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