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) SETERRQ2(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 = NULL; 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=NULL, pmat=NULL; 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 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetCoarseningComponent_C",NULL);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 243 PetscErrorCode PCView_HMG(PC pc,PetscViewer viewer) 244 { 245 PC_MG *mg = (PC_MG*)pc->data; 246 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 247 PetscErrorCode ierr; 248 PetscBool iascii; 249 250 PetscFunctionBegin; 251 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 252 if (iascii) { 253 ierr = PetscViewerASCIIPrintf(viewer," Reuse interpolation: %s\n",hmg->reuseinterp? "true":"false");CHKERRQ(ierr); 254 ierr = PetscViewerASCIIPrintf(viewer," Use subspace coarsening: %s\n",hmg->subcoarsening? "true":"false");CHKERRQ(ierr); 255 ierr = PetscViewerASCIIPrintf(viewer," Coarsening component: %D \n",hmg->component);CHKERRQ(ierr); 256 ierr = PetscViewerASCIIPrintf(viewer," Use MatMAIJ: %s \n",hmg->usematmaij? "true":"false");CHKERRQ(ierr); 257 ierr = PetscViewerASCIIPrintf(viewer," Inner PC type: %s \n",hmg->innerpctype);CHKERRQ(ierr); 258 } 259 ierr = PCView_MG(pc,viewer);CHKERRQ(ierr); 260 PetscFunctionReturn(0); 261 } 262 263 PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems *PetscOptionsObject,PC pc) 264 { 265 PetscErrorCode ierr; 266 PC_MG *mg = (PC_MG*)pc->data; 267 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 268 269 PetscFunctionBegin; 270 ierr = PetscOptionsHead(PetscOptionsObject,"HMG");CHKERRQ(ierr); 271 ierr = PetscOptionsBool("-pc_hmg_reuse_interpolation","Reuse the interpolation operators when possible (cheaper, weaker when matrix entries change a lot)","PCHMGSetReuseInterpolation",hmg->reuseinterp,&hmg->reuseinterp,NULL);CHKERRQ(ierr); 272 ierr = PetscOptionsBool("-pc_hmg_use_subspace_coarsening","Use the subspace coarsening to compute the interpolations","PCHMGSetUseSubspaceCoarsening",hmg->subcoarsening,&hmg->subcoarsening,NULL);CHKERRQ(ierr); 273 ierr = PetscOptionsBool("-pc_hmg_use_matmaij","Use MatMAIJ store interpolation for saving memory","PCHMGSetInnerPCType",hmg->usematmaij,&hmg->usematmaij,NULL);CHKERRQ(ierr); 274 ierr = PetscOptionsInt("-pc_hmg_coarsening_component","Which component is chosen for the subspace-based coarsening algorithm","PCHMGSetCoarseningComponent",hmg->component,&hmg->component,NULL);CHKERRQ(ierr); 275 ierr = PetscOptionsTail();CHKERRQ(ierr); 276 PetscFunctionReturn(0); 277 } 278 279 static PetscErrorCode PCHMGSetReuseInterpolation_HMG(PC pc, PetscBool reuse) 280 { 281 PC_MG *mg = (PC_MG*)pc->data; 282 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 283 284 PetscFunctionBegin; 285 hmg->reuseinterp = reuse; 286 PetscFunctionReturn(0); 287 } 288 289 /*@ 290 PCHMGSetReuseInterpolation - Reuse interpolation matrices in HMG 291 292 Logically Collective on PC 293 294 Input Parameters: 295 + pc - the HMG context 296 - reuse - True indicates that HMG will reuse the interpolations 297 298 Options Database Keys: 299 . -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time. 300 301 Level: beginner 302 303 .keywords: HMG, multigrid, interpolation, reuse, set 304 305 .seealso: PCHMG 306 @*/ 307 PetscErrorCode PCHMGSetReuseInterpolation(PC pc, PetscBool reuse) 308 { 309 PetscErrorCode ierr; 310 311 PetscFunctionBegin; 312 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 313 ierr = PetscUseMethod(pc,"PCHMGSetReuseInterpolation_C",(PC,PetscBool),(pc,reuse));CHKERRQ(ierr); 314 PetscFunctionReturn(0); 315 } 316 317 static PetscErrorCode PCHMGSetUseSubspaceCoarsening_HMG(PC pc, PetscBool subspace) 318 { 319 PC_MG *mg = (PC_MG*)pc->data; 320 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 321 322 PetscFunctionBegin; 323 hmg->subcoarsening = subspace; 324 PetscFunctionReturn(0); 325 } 326 327 /*@ 328 PCHMGSetUseSubspaceCoarsening - Use subspace coarsening in HMG 329 330 Logically Collective on PC 331 332 Input Parameters: 333 + pc - the HMG context 334 - reuse - True indicates that HMG will use the subspace coarsening 335 336 Options Database Keys: 337 . -pc_hmg_use_subspace_coarsening <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix). 338 339 Level: beginner 340 341 .keywords: HMG, multigrid, interpolation, subspace, coarsening 342 343 .seealso: PCHMG 344 @*/ 345 PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC pc, PetscBool subspace) 346 { 347 PetscErrorCode ierr; 348 349 PetscFunctionBegin; 350 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 351 ierr = PetscUseMethod(pc,"PCHMGSetUseSubspaceCoarsening_C",(PC,PetscBool),(pc,subspace));CHKERRQ(ierr); 352 PetscFunctionReturn(0); 353 } 354 355 static PetscErrorCode PCHMGSetInnerPCType_HMG(PC pc, PCType type) 356 { 357 PC_MG *mg = (PC_MG*)pc->data; 358 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 359 PetscErrorCode ierr; 360 361 PetscFunctionBegin; 362 ierr = PetscStrallocpy(type,&(hmg->innerpctype));CHKERRQ(ierr); 363 PetscFunctionReturn(0); 364 } 365 366 /*@C 367 PCHMGSetInnerPCType - Set an inner PC type 368 369 Logically Collective on PC 370 371 Input Parameters: 372 + pc - the HMG context 373 - type - <hypre, gamg> coarsening algorithm 374 375 Options Database Keys: 376 . -hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix 377 378 Level: beginner 379 380 .keywords: HMG, multigrid, interpolation, coarsening 381 382 .seealso: PCHMG, PCType 383 @*/ 384 PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type) 385 { 386 PetscErrorCode ierr; 387 388 PetscFunctionBegin; 389 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 390 ierr = PetscUseMethod(pc,"PCHMGSetInnerPCType_C",(PC,PCType),(pc,type));CHKERRQ(ierr); 391 PetscFunctionReturn(0); 392 } 393 394 static PetscErrorCode PCHMGSetCoarseningComponent_HMG(PC pc, PetscInt component) 395 { 396 PC_MG *mg = (PC_MG*)pc->data; 397 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 398 399 PetscFunctionBegin; 400 hmg->component = component; 401 PetscFunctionReturn(0); 402 } 403 404 /*@ 405 PCHMGSetCoarseningComponent - Set which component is used for the subspace-based coarsening algorithm 406 407 Logically Collective on PC 408 409 Input Parameters: 410 + pc - the HMG context 411 - component - which component PC will coarsen 412 413 Options Database Keys: 414 . -pc_hmg_coarsening_component - Which component is chosen for the subspace-based coarsening algorithm 415 416 Level: beginner 417 418 .keywords: HMG, multigrid, interpolation, coarsening, component 419 420 .seealso: PCHMG, PCType 421 @*/ 422 PetscErrorCode PCHMGSetCoarseningComponent(PC pc, PetscInt component) 423 { 424 PetscErrorCode ierr; 425 426 PetscFunctionBegin; 427 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 428 ierr = PetscUseMethod(pc,"PCHMGSetCoarseningComponent_C",(PC,PetscInt),(pc,component));CHKERRQ(ierr); 429 PetscFunctionReturn(0); 430 } 431 432 static PetscErrorCode PCHMGUseMatMAIJ_HMG(PC pc, PetscBool usematmaij) 433 { 434 PC_MG *mg = (PC_MG*)pc->data; 435 PC_HMG *hmg = (PC_HMG*) mg->innerctx; 436 437 PetscFunctionBegin; 438 hmg->usematmaij = usematmaij; 439 PetscFunctionReturn(0); 440 } 441 442 /*@ 443 PCHMGUseMatMAIJ - Set a flag that indicates if or not to use MatMAIJ for interpolations for saving memory 444 445 Logically Collective on PC 446 447 Input Parameters: 448 + pc - the HMG context 449 - usematmaij - if or not to use MatMAIJ for interpolations. By default, it is true for saving memory 450 451 Options Database Keys: 452 . -pc_hmg_use_matmaij - <true | false > 453 454 Level: beginner 455 456 .keywords: HMG, multigrid, interpolation, coarsening, MatMAIJ 457 458 .seealso: PCHMG, PCType 459 @*/ 460 PetscErrorCode PCHMGUseMatMAIJ(PC pc, PetscBool usematmaij) 461 { 462 PetscErrorCode ierr; 463 464 PetscFunctionBegin; 465 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 466 ierr = PetscUseMethod(pc,"PCHMGUseMatMAIJ_C",(PC,PetscBool),(pc,usematmaij));CHKERRQ(ierr); 467 PetscFunctionReturn(0); 468 } 469 470 /*MC 471 PCHMG - Hybrid of PETSc preconditioners (such as ASM, BJacobi, SOR, etc.) and Hypre BoomerAMG, GAMG or other multilevel methods. BoomerAMG, GAMG 472 or other multilevel methods is used to coarsen matrix and generate a sequence of coarse matrices and interpolations. The matrices and 473 interpolations are employed to construct PCMG, and then any available PETSc preconditioners can be chosen as smoothers and the coarse solver. 474 475 Options Database Keys: 476 + -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time. 477 . -pc_hmg_use_subspace_coarsening <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix). 478 . -hmg_inner_pc_type <hypre, gamg, ...> - What method is used to coarsen matrix 479 - -pc_hmg_use_matmaij <true | false> - Whether or not to use MatMAIJ for multicomponent problems for saving memory 480 481 482 Notes: 483 For multicomponent problems, we can just coarsen one submatrix associated with one particular component. In this way, the preconditioner setup 484 time is significantly reduced. One typical use case is neutron transport equations. There are many variables on each mesh vertex due to the 485 of angle and energy. Each variable, in fact, corresponds to the same PDEs but with different material properties. 486 487 Level: beginner 488 489 Concepts: Hybrid of ASM and MG, Subspace Coarsening 490 491 References: 492 . 1. - Fande Kong, Yaqi Wang, Derek R Gaston, Cody J Permann, Andrew E Slaughter, Alexander D Lindsay, Richard C Martineau, A highly parallel multilevel 493 Newton-Krylov-Schwarz method with subspace-based coarsening and partition-based balancing for the multigroup neutron transport equations on 494 3D unstructured meshes, arXiv preprint arXiv:1903.03659, 2019 495 496 .seealso: PCCreate(), PCSetType(), PCType, PC, PCMG, PCHYPRE, PCHMG, PCGetCoarseOperators(), PCGetInterpolations(), PCHMGSetReuseInterpolation(), PCHMGSetUseSubspaceCoarsening(), 497 PCHMGSetInnerPCType() 498 499 M*/ 500 PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc) 501 { 502 PetscErrorCode ierr; 503 PC_HMG *hmg; 504 PC_MG *mg; 505 506 PetscFunctionBegin; 507 /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */ 508 if (pc->ops->destroy) { 509 ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 510 pc->data = NULL; 511 } 512 ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr); 513 514 ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); 515 ierr = PetscObjectChangeTypeName((PetscObject)pc, PCHMG);CHKERRQ(ierr); 516 ierr = PetscNew(&hmg);CHKERRQ(ierr); 517 518 mg = (PC_MG*) pc->data; 519 mg->innerctx = hmg; 520 hmg->reuseinterp = PETSC_FALSE; 521 hmg->subcoarsening = PETSC_FALSE; 522 hmg->usematmaij = PETSC_TRUE; 523 hmg->component = 0; 524 hmg->innerpc = NULL; 525 526 pc->ops->setfromoptions = PCSetFromOptions_HMG; 527 pc->ops->view = PCView_HMG; 528 pc->ops->destroy = PCDestroy_HMG; 529 pc->ops->setup = PCSetUp_HMG; 530 531 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",PCHMGSetReuseInterpolation_HMG);CHKERRQ(ierr); 532 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",PCHMGSetUseSubspaceCoarsening_HMG);CHKERRQ(ierr); 533 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",PCHMGSetInnerPCType_HMG);CHKERRQ(ierr); 534 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetCoarseningComponent_C",PCHMGSetCoarseningComponent_HMG);CHKERRQ(ierr); 535 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGUseMatMAIJ_C",PCHMGUseMatMAIJ_HMG);CHKERRQ(ierr); 536 PetscFunctionReturn(0); 537 } 538