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