1360ee056SFande Kong #include <petscdm.h> 2360ee056SFande Kong #include <petscctable.h> 3360ee056SFande Kong #include <petsc/private/matimpl.h> 4360ee056SFande Kong #include <petsc/private/pcmgimpl.h> 5360ee056SFande Kong #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 6360ee056SFande Kong 7360ee056SFande Kong typedef struct { 88a2c336bSFande Kong PC innerpc; /* A MG inner PC (Hypre or PCGAMG) to setup interpolations and coarse operators */ 98a2c336bSFande Kong char* innerpctype; /* PCGAMG or PCHYPRE */ 108a2c336bSFande Kong PetscBool reuseinterp; /* A flag indicates if or not to reuse the interpolations */ 118a2c336bSFande Kong PetscBool subcoarsening; 124bb91820SFande Kong PetscBool usematmaij; 13360ee056SFande Kong } PC_HMG; 14360ee056SFande Kong 15360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems*,PC); 16360ee056SFande Kong PetscErrorCode PCReset_MG(PC); 178a2c336bSFande Kong 188a2c336bSFande Kong static PetscErrorCode PCHMGExtractSubMatrix_Private(Mat pmat,Mat *submat,MatReuse reuse,PetscInt blocksize) 19360ee056SFande Kong { 20360ee056SFande Kong IS isrow; 21360ee056SFande Kong PetscErrorCode ierr; 228a2c336bSFande Kong PetscInt rstart,rend; 23360ee056SFande Kong MPI_Comm comm; 24360ee056SFande Kong 25360ee056SFande Kong PetscFunctionBegin; 26360ee056SFande Kong ierr = PetscObjectGetComm((PetscObject)pmat,&comm);CHKERRQ(ierr); 27360ee056SFande Kong ierr = MatGetOwnershipRange(pmat,&rstart,&rend);CHKERRQ(ierr); 28360ee056SFande Kong if ((rend-rstart)%blocksize != 0) SETERRQ3(comm,PETSC_ERR_ARG_INCOMP,"Block size %d is inconsisent for [%d, %d) \n",blocksize,rstart,rend); 298a2c336bSFande Kong ierr = ISCreateStride(comm,(rend-rstart)/blocksize,rstart,blocksize,&isrow); 30360ee056SFande Kong ierr = MatCreateSubMatrix(pmat,isrow,isrow,reuse,submat);CHKERRQ(ierr); 31360ee056SFande Kong ierr = ISDestroy(&isrow);CHKERRQ(ierr); 32360ee056SFande Kong PetscFunctionReturn(0); 33360ee056SFande Kong } 34360ee056SFande Kong 358a2c336bSFande Kong static PetscErrorCode PCHMGExpandInterpolation_Private(Mat subinterp, Mat *interp, PetscInt blocksize) 36360ee056SFande Kong { 37360ee056SFande Kong PetscInt subrstart,subrend,subrowsize,subcolsize,subcstart,subcend,rowsize,colsize; 388a2c336bSFande Kong PetscInt subrow,row,nz,*d_nnz,*o_nnz,i,j,dnz,onz,max_nz,*indices; 398a2c336bSFande Kong const PetscInt *idx; 408a2c336bSFande Kong const PetscScalar *values; 41360ee056SFande Kong PetscErrorCode ierr; 428a2c336bSFande Kong MPI_Comm comm; 43360ee056SFande Kong 44360ee056SFande Kong PetscFunctionBegin; 458a2c336bSFande Kong ierr = PetscObjectGetComm((PetscObject)subinterp,&comm);CHKERRQ(ierr); 468a2c336bSFande Kong ierr = MatGetOwnershipRange(subinterp,&subrstart,&subrend);CHKERRQ(ierr); 478a2c336bSFande Kong subrowsize = subrend-subrstart; 48360ee056SFande Kong rowsize = subrowsize*blocksize; 49360ee056SFande Kong ierr = PetscCalloc2(rowsize,&d_nnz,rowsize,&o_nnz);CHKERRQ(ierr); 508a2c336bSFande Kong ierr = MatGetOwnershipRangeColumn(subinterp,&subcstart,&subcend);CHKERRQ(ierr); 518a2c336bSFande Kong subcolsize = subcend - subcstart; 528a2c336bSFande Kong colsize = subcolsize*blocksize; 53360ee056SFande Kong max_nz = 0; 548a2c336bSFande Kong for (subrow=subrstart;subrow<subrend;subrow++) { 558a2c336bSFande Kong ierr = MatGetRow(subinterp,subrow,&nz,&idx,NULL);CHKERRQ(ierr); 56360ee056SFande Kong if (max_nz<nz) max_nz = nz; 57360ee056SFande Kong dnz = 0; onz = 0; 58360ee056SFande Kong for (i=0;i<nz;i++) { 598a2c336bSFande Kong if(idx[i]>=subcstart && idx[i]<subcend) dnz++; 608a2c336bSFande Kong else onz++; 61360ee056SFande Kong } 62360ee056SFande Kong for (i=0;i<blocksize;i++) { 63360ee056SFande Kong d_nnz[(subrow-subrstart)*blocksize+i] = dnz; 64360ee056SFande Kong o_nnz[(subrow-subrstart)*blocksize+i] = onz; 65360ee056SFande Kong } 668a2c336bSFande Kong ierr = MatRestoreRow(subinterp,subrow,&nz,&idx,NULL);CHKERRQ(ierr); 67360ee056SFande Kong } 68360ee056SFande Kong ierr = MatCreateAIJ(comm,rowsize,colsize,PETSC_DETERMINE,PETSC_DETERMINE,0,d_nnz,0,o_nnz,interp);CHKERRQ(ierr); 69360ee056SFande Kong ierr = MatSetOption(*interp,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 70360ee056SFande Kong ierr = MatSetOption(*interp,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 71360ee056SFande Kong ierr = MatSetOption(*interp,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 72360ee056SFande Kong ierr = MatSetFromOptions(*interp);CHKERRQ(ierr); 73360ee056SFande Kong 74360ee056SFande Kong ierr = MatSetUp(*interp);CHKERRQ(ierr); 75360ee056SFande Kong ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 76071fcb05SBarry Smith ierr = PetscMalloc1(max_nz,&indices);CHKERRQ(ierr); 778a2c336bSFande Kong for (subrow=subrstart; subrow<subrend; subrow++) { 788a2c336bSFande Kong ierr = MatGetRow(subinterp,subrow,&nz,&idx,&values);CHKERRQ(ierr); 79360ee056SFande Kong for (i=0;i<blocksize;i++) { 80360ee056SFande Kong row = subrow*blocksize+i; 81360ee056SFande Kong for (j=0;j<nz;j++) { 82360ee056SFande Kong indices[j] = idx[j]*blocksize+i; 83360ee056SFande Kong } 84360ee056SFande Kong ierr = MatSetValues(*interp,1,&row,nz,indices,values,INSERT_VALUES);CHKERRQ(ierr); 85360ee056SFande Kong } 868a2c336bSFande Kong ierr = MatRestoreRow(subinterp,subrow,&nz,&idx,&values);CHKERRQ(ierr); 87360ee056SFande Kong } 88360ee056SFande Kong ierr = PetscFree(indices);CHKERRQ(ierr); 89360ee056SFande Kong ierr = MatAssemblyBegin(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 90360ee056SFande Kong ierr = MatAssemblyEnd(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 91360ee056SFande Kong PetscFunctionReturn(0); 92360ee056SFande Kong } 93360ee056SFande Kong 94360ee056SFande Kong PetscErrorCode PCSetUp_HMG(PC pc) 95360ee056SFande Kong { 96360ee056SFande Kong PetscErrorCode ierr; 97360ee056SFande Kong Mat PA, submat; 98360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 99360ee056SFande Kong PC_HMG *hmg = (PC_HMG*) mg->innerctx; 100360ee056SFande Kong MPI_Comm comm; 101360ee056SFande Kong PetscInt level; 102360ee056SFande Kong PetscInt num_levels; 1038a2c336bSFande Kong Mat *operators,*interpolations; 1048a2c336bSFande Kong PetscInt blocksize; 105fd2dd295SFande Kong const char *prefix; 10607a4832bSFande Kong PCMGGalerkinType galerkin; 107360ee056SFande Kong 108360ee056SFande Kong PetscFunctionBegin; 109360ee056SFande Kong ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 110360ee056SFande Kong if (pc->setupcalled) { 1115f36d8f2SFande Kong if (hmg->reuseinterp) { 1125f36d8f2SFande Kong /* If we did not use Galerkin in the last call or we have a different sparsity pattern now, 1135f36d8f2SFande Kong * we have to build from scratch 11407a4832bSFande Kong * */ 11507a4832bSFande Kong ierr = PCMGGetGalerkin(pc,&galerkin);CHKERRQ(ierr); 1165f36d8f2SFande Kong if (galerkin == PC_MG_GALERKIN_NONE || pc->flag != SAME_NONZERO_PATTERN) pc->setupcalled = PETSC_FALSE; 117360ee056SFande Kong ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT);CHKERRQ(ierr); 118360ee056SFande Kong ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 119360ee056SFande Kong PetscFunctionReturn(0); 120360ee056SFande Kong } else { 121360ee056SFande Kong ierr = PCReset_MG(pc);CHKERRQ(ierr); 122360ee056SFande Kong pc->setupcalled = PETSC_FALSE; 123360ee056SFande Kong } 124360ee056SFande Kong } 125360ee056SFande Kong 1268a2c336bSFande Kong /* Create an inner PC (GAMG or HYPRE) */ 1278a2c336bSFande Kong if (!hmg->innerpc) { 1288a2c336bSFande Kong ierr = PCCreate(comm,&hmg->innerpc);CHKERRQ(ierr); 129fd2dd295SFande Kong /* If users do not set an inner pc type, we need to set a default value */ 130fd2dd295SFande Kong if (!hmg->innerpctype) { 131fd2dd295SFande Kong /* If hypre is available, use hypre, otherwise, use gamg */ 132fd2dd295SFande Kong #if PETSC_HAVE_HYPRE 133fd2dd295SFande Kong ierr = PetscStrallocpy(PCHYPRE,&(hmg->innerpctype));CHKERRQ(ierr); 134fd2dd295SFande Kong #else 135fd2dd295SFande Kong ierr = PetscStrallocpy(PCGAMG,&(hmg->innerpctype));CHKERRQ(ierr); 136fd2dd295SFande Kong #endif 1378a2c336bSFande Kong } 138fd2dd295SFande Kong ierr = PCSetType(hmg->innerpc,hmg->innerpctype);CHKERRQ(ierr); 139360ee056SFande Kong } 140360ee056SFande Kong ierr = PCGetOperators(pc,NULL,&PA);CHKERRQ(ierr); 1418a2c336bSFande Kong /* Users need to correctly set a block size of matrix in order to use subspace coarsening */ 1428a2c336bSFande Kong ierr = MatGetBlockSize(PA,&blocksize);CHKERRQ(ierr); 1438a2c336bSFande Kong if (blocksize<=1) hmg->subcoarsening = PETSC_FALSE; 1448a2c336bSFande Kong /* Extract a submatrix for constructing subinterpolations */ 1458a2c336bSFande Kong if (hmg->subcoarsening) { 1468a2c336bSFande Kong ierr = PCHMGExtractSubMatrix_Private(PA,&submat,MAT_INITIAL_MATRIX,blocksize);CHKERRQ(ierr); 147360ee056SFande Kong PA = submat; 148360ee056SFande Kong } 1498a2c336bSFande Kong ierr = PCSetOperators(hmg->innerpc,PA,PA);CHKERRQ(ierr); 1508a2c336bSFande Kong if (hmg->subcoarsening) { 151360ee056SFande Kong ierr = MatDestroy(&PA);CHKERRQ(ierr); 152360ee056SFande Kong } 1538a2c336bSFande Kong /* Setup inner PC correctly. During this step, matrix will be coarsened */ 1548a2c336bSFande Kong ierr = PCSetUseAmat(hmg->innerpc,PETSC_FALSE);CHKERRQ(ierr); 155fd2dd295SFande Kong ierr = PetscObjectGetOptionsPrefix((PetscObject)pc,&prefix);CHKERRQ(ierr); 156fd2dd295SFande Kong ierr = PetscObjectSetOptionsPrefix((PetscObject)hmg->innerpc,prefix);CHKERRQ(ierr); 157fd2dd295SFande Kong ierr = PetscObjectAppendOptionsPrefix((PetscObject)hmg->innerpc,"hmg_inner_");CHKERRQ(ierr); 158fd2dd295SFande Kong ierr = PCSetFromOptions(hmg->innerpc);CHKERRQ(ierr); 1598a2c336bSFande Kong ierr = PCSetUp(hmg->innerpc); 160360ee056SFande Kong 1618a2c336bSFande Kong /* Obtain interpolations IN PLACE. For BoomerAMG, (I,J,data) is reused to avoid memory overhead */ 162fd2dd295SFande Kong ierr = PCGetInterpolations(hmg->innerpc,&num_levels,&interpolations);CHKERRQ(ierr); 1638a2c336bSFande Kong /* We can reuse the coarse operators when we do the full space coarsening */ 1648a2c336bSFande Kong if (!hmg->subcoarsening) { 165fd2dd295SFande Kong ierr = PCGetCoarseOperators(hmg->innerpc,&num_levels,&operators);CHKERRQ(ierr); 166360ee056SFande Kong } 167360ee056SFande Kong 1688a2c336bSFande Kong ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr); 1698a2c336bSFande Kong hmg->innerpc = 0; 1708a2c336bSFande Kong ierr = PCMGSetLevels_MG(pc,num_levels,NULL);CHKERRQ(ierr); 1718a2c336bSFande Kong /* Set coarse matrices and interpolations to PCMG */ 172360ee056SFande Kong for (level=num_levels-1; level>0; level--) { 173360ee056SFande Kong Mat P=0, pmat=0; 174360ee056SFande Kong Vec b, x,r; 1758a2c336bSFande Kong if (hmg->subcoarsening) { 1764bb91820SFande Kong if (hmg->usematmaij) { 1774bb91820SFande Kong ierr = MatCreateMAIJ(interpolations[level-1],blocksize,&P);CHKERRQ(ierr); 1784bb91820SFande Kong ierr = MatDestroy(&interpolations[level-1]);CHKERRQ(ierr); 1794bb91820SFande Kong } else { 1808a2c336bSFande Kong /* Grow interpolation. In the future, we should use MAIJ */ 1818a2c336bSFande Kong ierr = PCHMGExpandInterpolation_Private(interpolations[level-1],&P,blocksize);CHKERRQ(ierr); 1828a2c336bSFande Kong ierr = MatDestroy(&interpolations[level-1]);CHKERRQ(ierr); 1834bb91820SFande Kong } 184360ee056SFande Kong } else { 1858a2c336bSFande Kong P = interpolations[level-1]; 186360ee056SFande Kong } 187360ee056SFande Kong ierr = MatCreateVecs(P,&b,&r);CHKERRQ(ierr); 188360ee056SFande Kong ierr = PCMGSetInterpolation(pc,level,P);CHKERRQ(ierr); 189360ee056SFande Kong ierr = PCMGSetRestriction(pc,level,P);CHKERRQ(ierr); 190360ee056SFande Kong ierr = MatDestroy(&P);CHKERRQ(ierr); 1918a2c336bSFande Kong /* We reuse the matrices when we do not do subspace coarsening */ 1928a2c336bSFande Kong if ((level-1)>=0 && !hmg->subcoarsening) { 1938a2c336bSFande Kong pmat = operators[level-1]; 1948a2c336bSFande Kong ierr = PCMGSetOperators(pc,level-1,pmat,pmat);CHKERRQ(ierr); 195360ee056SFande Kong ierr = MatDestroy(&pmat);CHKERRQ(ierr); 196360ee056SFande Kong } 197360ee056SFande Kong ierr = PCMGSetRhs(pc,level-1,b);CHKERRQ(ierr); 198360ee056SFande Kong 199360ee056SFande Kong ierr = PCMGSetR(pc,level,r);CHKERRQ(ierr); 200360ee056SFande Kong ierr = VecDestroy(&r);CHKERRQ(ierr); 201360ee056SFande Kong 202fd2dd295SFande Kong ierr = VecDuplicate(b,&x);CHKERRQ(ierr); 203360ee056SFande Kong ierr = PCMGSetX(pc,level-1,x);CHKERRQ(ierr); 204360ee056SFande Kong ierr = VecDestroy(&x);CHKERRQ(ierr); 205360ee056SFande Kong ierr = VecDestroy(&b);CHKERRQ(ierr); 206360ee056SFande Kong } 2078a2c336bSFande Kong ierr = PetscFree(interpolations);CHKERRQ(ierr); 2088a2c336bSFande Kong if (!hmg->subcoarsening) { 2098a2c336bSFande Kong ierr = PetscFree(operators);CHKERRQ(ierr); 2108a2c336bSFande Kong } 2118a2c336bSFande Kong /* Turn Galerkin off when we already have coarse operators */ 2128a2c336bSFande Kong ierr = PCMGSetGalerkin(pc,hmg->subcoarsening ? PC_MG_GALERKIN_PMAT:PC_MG_GALERKIN_NONE);CHKERRQ(ierr); 213360ee056SFande Kong ierr = PCSetDM(pc,NULL);CHKERRQ(ierr); 214360ee056SFande Kong ierr = PCSetUseAmat(pc,PETSC_FALSE);CHKERRQ(ierr); 215360ee056SFande Kong ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 216360ee056SFande Kong ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */ 217360ee056SFande Kong ierr = PetscOptionsEnd();CHKERRQ(ierr); 218360ee056SFande Kong ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 219360ee056SFande Kong PetscFunctionReturn(0); 220360ee056SFande Kong } 221360ee056SFande Kong 222360ee056SFande Kong PetscErrorCode PCDestroy_HMG(PC pc) 223360ee056SFande Kong { 224360ee056SFande Kong PetscErrorCode ierr; 225360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 2268a2c336bSFande Kong PC_HMG *hmg = (PC_HMG*) mg->innerctx; 227360ee056SFande Kong 228360ee056SFande Kong PetscFunctionBegin; 2298a2c336bSFande Kong ierr = PCDestroy(&hmg->innerpc);CHKERRQ(ierr); 2308a2c336bSFande Kong ierr = PetscFree(hmg->innerpctype);CHKERRQ(ierr); 2318a2c336bSFande Kong ierr = PetscFree(hmg);CHKERRQ(ierr); 232360ee056SFande Kong ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 233fd2dd295SFande Kong 234fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",NULL);CHKERRQ(ierr); 235fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",NULL);CHKERRQ(ierr); 236fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",NULL);CHKERRQ(ierr); 237360ee056SFande Kong PetscFunctionReturn(0); 238360ee056SFande Kong } 239360ee056SFande Kong 240360ee056SFande Kong PetscErrorCode PCView_HMG(PC pc,PetscViewer viewer) 241360ee056SFande Kong { 242360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 2438a2c336bSFande Kong PC_HMG *hmg = (PC_HMG*) mg->innerctx; 244360ee056SFande Kong PetscErrorCode ierr; 245360ee056SFande Kong PetscBool iascii; 246360ee056SFande Kong 247360ee056SFande Kong PetscFunctionBegin; 248360ee056SFande Kong ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 249360ee056SFande Kong if (iascii) { 2508a2c336bSFande Kong ierr = PetscViewerASCIIPrintf(viewer," Reuse interpolation: %s\n",hmg->reuseinterp? "true":"false");CHKERRQ(ierr); 2518a2c336bSFande Kong ierr = PetscViewerASCIIPrintf(viewer," Use subspace coarsening: %s\n",hmg->subcoarsening? "true":"false");CHKERRQ(ierr); 2528a2c336bSFande Kong ierr = PetscViewerASCIIPrintf(viewer," Inner PC type: %s \n",hmg->innerpctype);CHKERRQ(ierr); 253360ee056SFande Kong } 254360ee056SFande Kong ierr = PCView_MG(pc,viewer);CHKERRQ(ierr); 255360ee056SFande Kong PetscFunctionReturn(0); 256360ee056SFande Kong } 257360ee056SFande Kong 258360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems *PetscOptionsObject,PC pc) 259360ee056SFande Kong { 260360ee056SFande Kong PetscErrorCode ierr; 261360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 2628a2c336bSFande Kong PC_HMG *hmg = (PC_HMG*) mg->innerctx; 263360ee056SFande Kong 264360ee056SFande Kong PetscFunctionBegin; 265360ee056SFande Kong ierr = PetscOptionsHead(PetscOptionsObject,"HMG");CHKERRQ(ierr); 2668a2c336bSFande Kong 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); 2678a2c336bSFande Kong ierr = PetscOptionsBool("-pc_hmg_use_subspace_coarsening","Use the subspace coarsening to compute the interpolations","None",hmg->subcoarsening,&hmg->subcoarsening,NULL);CHKERRQ(ierr); 2684bb91820SFande Kong ierr = PetscOptionsBool("-pc_hmg_use_matmaij","Use MatMAIJ store interpolation for saving memory","None",hmg->usematmaij,&hmg->usematmaij,NULL);CHKERRQ(ierr); 269360ee056SFande Kong ierr = PetscOptionsTail();CHKERRQ(ierr); 270360ee056SFande Kong PetscFunctionReturn(0); 271360ee056SFande Kong } 272360ee056SFande Kong 273fd2dd295SFande Kong static PetscErrorCode PCHMGSetReuseInterpolation_HMG(PC pc, PetscBool reuse) 274fd2dd295SFande Kong { 27507a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 27607a4832bSFande Kong PC_HMG *hmg = (PC_HMG*) mg->innerctx; 277fd2dd295SFande Kong 278fd2dd295SFande Kong PetscFunctionBegin; 279fd2dd295SFande Kong hmg->reuseinterp = reuse; 280fd2dd295SFande Kong PetscFunctionReturn(0); 281fd2dd295SFande Kong } 282fd2dd295SFande Kong 283360ee056SFande Kong /*MC 2848a2c336bSFande Kong PCHMGSetReuseInterpolation - Reuse interpolation matrices in HMG 2858a2c336bSFande Kong 2868a2c336bSFande Kong Logically Collective on PC 2878a2c336bSFande Kong 2888a2c336bSFande Kong Input Parameters: 2898a2c336bSFande Kong + pc - the HMG context 2908a2c336bSFande Kong - reuse - True indicates that HMG will reuse the interpolations 2918a2c336bSFande Kong 2928a2c336bSFande Kong Options Database Keys: 2938a2c336bSFande Kong + -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time. 2948a2c336bSFande Kong 2958a2c336bSFande Kong Level: beginner 2968a2c336bSFande Kong 2978a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, reuse, set 2988a2c336bSFande Kong 2998a2c336bSFande Kong .seealso: PCHMG 3008a2c336bSFande Kong M*/ 3018a2c336bSFande Kong PetscErrorCode PCHMGSetReuseInterpolation(PC pc, PetscBool reuse) 3028a2c336bSFande Kong { 303fd2dd295SFande Kong PetscErrorCode ierr; 304fd2dd295SFande Kong 305fd2dd295SFande Kong PetscFunctionBegin; 306fd2dd295SFande Kong PetscValidHeaderSpecific(pc,PC_CLASSID,1); 307fd2dd295SFande Kong ierr = PetscUseMethod(pc,"PCHMGSetReuseInterpolation_C",(PC,PetscBool),(pc,reuse));CHKERRQ(ierr); 308fd2dd295SFande Kong PetscFunctionReturn(0); 309fd2dd295SFande Kong } 310fd2dd295SFande Kong 311fd2dd295SFande Kong static PetscErrorCode PCHMGSetUseSubspaceCoarsening_HMG(PC pc, PetscBool subspace) 312fd2dd295SFande Kong { 31307a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 31407a4832bSFande Kong PC_HMG *hmg = (PC_HMG*) mg->innerctx; 3158a2c336bSFande Kong 3168a2c336bSFande Kong PetscFunctionBegin; 317fd2dd295SFande Kong hmg->subcoarsening = subspace; 3188a2c336bSFande Kong PetscFunctionReturn(0); 3198a2c336bSFande Kong } 3208a2c336bSFande Kong 3218a2c336bSFande Kong /*MC 3228a2c336bSFande Kong PCHMGSetUseSubspaceCoarsening - Use subspace coarsening in HMG 3238a2c336bSFande Kong 3248a2c336bSFande Kong Logically Collective on PC 3258a2c336bSFande Kong 3268a2c336bSFande Kong Input Parameters: 3278a2c336bSFande Kong + pc - the HMG context 3288a2c336bSFande Kong - reuse - True indicates that HMG will use the subspace coarsening 3298a2c336bSFande Kong 3308a2c336bSFande Kong Options Database Keys: 3318a2c336bSFande Kong + -pc_hmg_use_subspace_coarsening <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix). 3328a2c336bSFande Kong 3338a2c336bSFande Kong Level: beginner 3348a2c336bSFande Kong 3358a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, subspace, coarsening 3368a2c336bSFande Kong 3378a2c336bSFande Kong .seealso: PCHMG 3388a2c336bSFande Kong M*/ 3398a2c336bSFande Kong PetscErrorCode PCHMGSetUseSubspaceCoarsening(PC pc, PetscBool subspace) 3408a2c336bSFande Kong { 341fd2dd295SFande Kong PetscErrorCode ierr; 3428a2c336bSFande Kong 3438a2c336bSFande Kong PetscFunctionBegin; 3448a2c336bSFande Kong PetscValidHeaderSpecific(pc,PC_CLASSID,1); 345fd2dd295SFande Kong ierr = PetscUseMethod(pc,"PCHMGSetUseSubspaceCoarsening_C",(PC,PetscBool),(pc,subspace));CHKERRQ(ierr); 346fd2dd295SFande Kong PetscFunctionReturn(0); 347fd2dd295SFande Kong } 348fd2dd295SFande Kong 349fd2dd295SFande Kong static PetscErrorCode PCHMGSetInnerPCType_HMG(PC pc, PCType type) 350fd2dd295SFande Kong { 35107a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 35207a4832bSFande Kong PC_HMG *hmg = (PC_HMG*) mg->innerctx; 353fd2dd295SFande Kong PetscErrorCode ierr; 354fd2dd295SFande Kong 355fd2dd295SFande Kong PetscFunctionBegin; 356fd2dd295SFande Kong ierr = PetscStrallocpy(type,&(hmg->innerpctype));CHKERRQ(ierr); 3578a2c336bSFande Kong PetscFunctionReturn(0); 3588a2c336bSFande Kong } 3598a2c336bSFande Kong 3608a2c336bSFande Kong /*MC 3618a2c336bSFande Kong PCHMGSetInnerPCType - Set an inner PC type 3628a2c336bSFande Kong 3638a2c336bSFande Kong Logically Collective on PC 3648a2c336bSFande Kong 3658a2c336bSFande Kong Input Parameters: 3668a2c336bSFande Kong + pc - the HMG context 3678a2c336bSFande Kong - type - <hypre, gamg> coarsening algorithm 3688a2c336bSFande Kong 3698a2c336bSFande Kong Options Database Keys: 370fd2dd295SFande Kong + -hmg_inner_pc_type <hypre, gamg> - What method is used to coarsen matrix 3718a2c336bSFande Kong 3728a2c336bSFande Kong Level: beginner 3738a2c336bSFande Kong 3748a2c336bSFande Kong .keywords: HMG, multigrid, interpolation, coarsening 3758a2c336bSFande Kong 3768a2c336bSFande Kong .seealso: PCHMG, PCType 3778a2c336bSFande Kong M*/ 3788a2c336bSFande Kong PetscErrorCode PCHMGSetInnerPCType(PC pc, PCType type) 3798a2c336bSFande Kong { 3808a2c336bSFande Kong PetscErrorCode ierr; 3818a2c336bSFande Kong 3828a2c336bSFande Kong PetscFunctionBegin; 3838a2c336bSFande Kong PetscValidHeaderSpecific(pc,PC_CLASSID,1); 384fd2dd295SFande Kong ierr = PetscUseMethod(pc,"PCHMGSetInnerPCType_C",(PC,PCType),(pc,type));CHKERRQ(ierr); 3858a2c336bSFande Kong PetscFunctionReturn(0); 3868a2c336bSFande Kong } 3878a2c336bSFande Kong 3888a2c336bSFande Kong /*MC 389fd2dd295SFande Kong PCHMG - Hybrid of PETSc preconditioners (such as ASM, BJacobi, SOR, etc.) and Hypre BoomerAMG, GAMG or other multilevel methods. BoomerAMG, GAMG 390fd2dd295SFande Kong or other multilevel methods is used to coarsen matrix and generate a sequence of coarse matrices and interpolations. The matrices and 391fd2dd295SFande Kong interpolations are employed to construct PCMG, and then any available PETSc preconditioners can be chosen as smoothers and the coarse solver. 392360ee056SFande Kong 393360ee056SFande Kong Options Database Keys: 3948a2c336bSFande Kong + -pc_hmg_reuse_interpolation <true | false> - Whether or not to reuse the interpolations. If true, it potentially save the compute time. 3958a2c336bSFande Kong . -pc_hmg_use_subspace_coarsening <true | false> - Whether or not to use subspace coarsening (that is, coarsen a submatrix). 396fd2dd295SFande Kong . -hmg_inner_pc_type <hypre, gamg, ...> - What method is used to coarsen matrix 3974bb91820SFande Kong - -pc_hmg_use_matmaij <true | false> - Whether or not to use MatMAIJ for multicomponent problems for saving memory 398360ee056SFande Kong 399360ee056SFande Kong 400360ee056SFande Kong Notes: 401360ee056SFande Kong For multicomponent problems, we can just coarsen one submatrix associated with one particular component. In this way, the preconditioner setup 402360ee056SFande Kong time is significantly reduced. One typical use case is neutron transport equations. There are many variables on each mesh vertex due to the 403360ee056SFande Kong of angle and energy. Each variable, in fact, corresponds to the same PDEs but with different material properties. 404360ee056SFande Kong 405360ee056SFande Kong Level: beginner 406360ee056SFande Kong 4078a2c336bSFande Kong Concepts: Hybrid of ASM and MG, Subspace Coarsening 408360ee056SFande Kong 409360ee056SFande Kong References: 410360ee056SFande Kong + 1. - Fande Kong, Yaqi Wang, Derek R Gaston, Cody J Permann, Andrew E Slaughter, Alexander D Lindsay, Richard C Martineau, A highly parallel multilevel 411360ee056SFande Kong Newton-Krylov-Schwarz method with subspace-based coarsening and partition-based balancing for the multigroup neutron transport equations on 412360ee056SFande Kong 3D unstructured meshes, arXiv preprint arXiv:1903.03659, 2019 413360ee056SFande Kong 414fd2dd295SFande Kong .seealso: PCCreate(), PCSetType(), PCType, PC, PCMG, PCHYPRE, PCHMG, PCGetCoarseOperators(), PCGetInterpolations(), PCHMGSetReuseInterpolation(), PCHMGSetUseSubspaceCoarsening(), 415fd2dd295SFande Kong PCHMGSetInnerPCType() 416360ee056SFande Kong 417360ee056SFande Kong M*/ 418360ee056SFande Kong PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc) 419360ee056SFande Kong { 420360ee056SFande Kong PetscErrorCode ierr; 421360ee056SFande Kong PC_HMG *hmg; 422360ee056SFande Kong PC_MG *mg; 423360ee056SFande Kong 424360ee056SFande Kong PetscFunctionBegin; 425360ee056SFande Kong /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */ 426360ee056SFande Kong if (pc->ops->destroy) { 427360ee056SFande Kong ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 428360ee056SFande Kong pc->data = 0; 429360ee056SFande Kong } 430360ee056SFande Kong ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr); 431360ee056SFande Kong 432360ee056SFande Kong ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); 433*645b8336SFande Kong ierr = PetscObjectChangeTypeName((PetscObject)pc, PCHMG);CHKERRQ(ierr); 434fd2dd295SFande Kong ierr = PetscNew(&hmg);CHKERRQ(ierr); 435360ee056SFande Kong 436360ee056SFande Kong mg = (PC_MG*) pc->data; 437360ee056SFande Kong mg->innerctx = hmg; 438360ee056SFande Kong hmg->reuseinterp = PETSC_FALSE; 4398a2c336bSFande Kong hmg->subcoarsening = PETSC_FALSE; 4404bb91820SFande Kong hmg->usematmaij = PETSC_TRUE; 4418a2c336bSFande Kong hmg->innerpc = NULL; 442360ee056SFande Kong 443360ee056SFande Kong pc->ops->setfromoptions = PCSetFromOptions_HMG; 444360ee056SFande Kong pc->ops->view = PCView_HMG; 445360ee056SFande Kong pc->ops->destroy = PCDestroy_HMG; 446360ee056SFande Kong pc->ops->setup = PCSetUp_HMG; 447fd2dd295SFande Kong 448fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetReuseInterpolation_C",PCHMGSetReuseInterpolation_HMG);CHKERRQ(ierr); 449fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetUseSubspaceCoarsening_C",PCHMGSetUseSubspaceCoarsening_HMG);CHKERRQ(ierr); 450fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCHMGSetInnerPCType_C",PCHMGSetInnerPCType_HMG);CHKERRQ(ierr); 451360ee056SFande Kong PetscFunctionReturn(0); 452360ee056SFande Kong } 453