1*360ee056SFande Kong 2*360ee056SFande Kong #include <petscdm.h> 3*360ee056SFande Kong #include <petscctable.h> 4*360ee056SFande Kong #include <petsc/private/matimpl.h> 5*360ee056SFande Kong /* Need to access the hypre private data */ 6*360ee056SFande Kong #include <_hypre_parcsr_ls.h> 7*360ee056SFande Kong #include <petsc/private/pcmgimpl.h> 8*360ee056SFande Kong #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 9*360ee056SFande Kong 10*360ee056SFande Kong typedef struct { 11*360ee056SFande Kong PC hypre; 12*360ee056SFande Kong PetscBool reuseinterp; 13*360ee056SFande Kong PetscInt blocksize; 14*360ee056SFande Kong } PC_HMG; 15*360ee056SFande Kong 16*360ee056SFande Kong PetscErrorCode MatConvert_ParCSRMatrix_AIJ(MPI_Comm, hypre_ParCSRMatrix*, MatType, MatReuse, Mat*); 17*360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems*,PC); 18*360ee056SFande Kong PetscErrorCode PCSetFromOptions_HYPRE(PetscOptionItems*,PC); 19*360ee056SFande Kong PetscErrorCode PCReset_MG(PC); 20*360ee056SFande Kong PetscErrorCode PCHYPREGetSolver(PC,HYPRE_Solver*); 21*360ee056SFande Kong 22*360ee056SFande Kong static PetscErrorCode PCHMGExtractSubMatrix_HMG(Mat pmat,Mat *submat,MatReuse reuse,PetscInt blocksize) 23*360ee056SFande Kong { 24*360ee056SFande Kong IS isrow; 25*360ee056SFande Kong PetscErrorCode ierr; 26*360ee056SFande Kong PetscInt rstart,rend, row,subsize, *rowsindices; 27*360ee056SFande Kong MPI_Comm comm; 28*360ee056SFande Kong 29*360ee056SFande Kong PetscFunctionBegin; 30*360ee056SFande Kong ierr = PetscObjectGetComm((PetscObject)pmat,&comm);CHKERRQ(ierr); 31*360ee056SFande Kong ierr = MatGetOwnershipRange(pmat,&rstart,&rend);CHKERRQ(ierr); 32*360ee056SFande Kong if ((rend-rstart)%blocksize != 0) SETERRQ3(comm,PETSC_ERR_ARG_INCOMP,"Block size %d is inconsisent for [%d, %d) \n",blocksize,rstart,rend); 33*360ee056SFande Kong subsize = (rend-rstart)/blocksize; 34*360ee056SFande Kong ierr = PetscCalloc1(subsize,&rowsindices);CHKERRQ(ierr); 35*360ee056SFande Kong subsize = 0; 36*360ee056SFande Kong for (row=rstart; row<rend; row+=blocksize){ 37*360ee056SFande Kong rowsindices[subsize++]=row; 38*360ee056SFande Kong } 39*360ee056SFande Kong ierr = ISCreateGeneral(comm,subsize,rowsindices,PETSC_OWN_POINTER,&isrow);CHKERRQ(ierr); 40*360ee056SFande Kong ierr = MatCreateSubMatrix(pmat,isrow,isrow,reuse,submat);CHKERRQ(ierr); 41*360ee056SFande Kong ierr = ISDestroy(&isrow);CHKERRQ(ierr); 42*360ee056SFande Kong PetscFunctionReturn(0); 43*360ee056SFande Kong } 44*360ee056SFande Kong 45*360ee056SFande Kong static PetscErrorCode PCHMGExpandInterpolation_HMG(MPI_Comm comm,hypre_ParCSRMatrix *subinterp, Mat *interp, PetscInt blocksize) 46*360ee056SFande Kong { 47*360ee056SFande Kong PetscInt subrstart,subrend,subrowsize,subcolsize,subcstart,subcend,rowsize,colsize; 48*360ee056SFande Kong PetscInt subrow,row,*idx,nz,*d_nnz,*o_nnz,i,j,dnz,onz,max_nz,*indices; 49*360ee056SFande Kong PetscScalar *values; 50*360ee056SFande Kong PetscErrorCode ierr; 51*360ee056SFande Kong 52*360ee056SFande Kong PetscFunctionBegin; 53*360ee056SFande Kong subrstart = hypre_ParCSRMatrixFirstRowIndex(subinterp); 54*360ee056SFande Kong subrend = hypre_ParCSRMatrixLastRowIndex(subinterp); 55*360ee056SFande Kong subrowsize = subrend-subrstart+1; 56*360ee056SFande Kong rowsize = subrowsize*blocksize; 57*360ee056SFande Kong subcstart = hypre_ParCSRMatrixFirstColDiag(subinterp); 58*360ee056SFande Kong subcend = hypre_ParCSRMatrixLastColDiag(subinterp); 59*360ee056SFande Kong subcolsize = subcend-subcstart+1; 60*360ee056SFande Kong colsize = subcolsize*blocksize; 61*360ee056SFande Kong ierr = PetscCalloc2(rowsize,&d_nnz,rowsize,&o_nnz);CHKERRQ(ierr); 62*360ee056SFande Kong max_nz = 0; 63*360ee056SFande Kong for(subrow=subrstart; subrow<=subrend; subrow++){ 64*360ee056SFande Kong HYPRE_ParCSRMatrixGetRow(subinterp,subrow,(HYPRE_Int*)&nz,(HYPRE_Int**)&idx,NULL); 65*360ee056SFande Kong if (max_nz<nz) max_nz = nz; 66*360ee056SFande Kong dnz = 0; onz = 0; 67*360ee056SFande Kong for(i=0;i<nz;i++){ 68*360ee056SFande Kong if(idx[i]<subcstart || idx[i]>subcend) onz++; 69*360ee056SFande Kong else dnz++; 70*360ee056SFande Kong } 71*360ee056SFande Kong for(i=0;i<blocksize;i++){ 72*360ee056SFande Kong d_nnz[(subrow-subrstart)*blocksize+i] = dnz; 73*360ee056SFande Kong o_nnz[(subrow-subrstart)*blocksize+i] = onz; 74*360ee056SFande Kong } 75*360ee056SFande Kong HYPRE_ParCSRMatrixRestoreRow(subinterp,subrow,(HYPRE_Int*)&nz,(HYPRE_Int**)&idx,NULL); 76*360ee056SFande Kong } 77*360ee056SFande Kong ierr = MatCreateAIJ(comm,rowsize,colsize,PETSC_DETERMINE,PETSC_DETERMINE,0,d_nnz,0,o_nnz,interp);CHKERRQ(ierr); 78*360ee056SFande Kong ierr = MatSetOption(*interp,MAT_IGNORE_OFF_PROC_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 79*360ee056SFande Kong ierr = MatSetOption(*interp,MAT_IGNORE_ZERO_ENTRIES,PETSC_TRUE);CHKERRQ(ierr); 80*360ee056SFande Kong ierr = MatSetOption(*interp,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 81*360ee056SFande Kong ierr = MatSetFromOptions(*interp);CHKERRQ(ierr); 82*360ee056SFande Kong 83*360ee056SFande Kong ierr = MatSetUp(*interp);CHKERRQ(ierr); 84*360ee056SFande Kong ierr = PetscFree2(d_nnz,o_nnz);CHKERRQ(ierr); 85*360ee056SFande Kong ierr = PetscCalloc1(max_nz,&indices);CHKERRQ(ierr); 86*360ee056SFande Kong for(subrow=subrstart; subrow<=subrend; subrow++){ 87*360ee056SFande Kong HYPRE_ParCSRMatrixGetRow(subinterp,subrow,(HYPRE_Int*)&nz,(HYPRE_Int**)&idx,&values); 88*360ee056SFande Kong for(i=0;i<blocksize;i++){ 89*360ee056SFande Kong row = subrow*blocksize+i; 90*360ee056SFande Kong for (j=0;j<nz;j++){ 91*360ee056SFande Kong indices[j] = idx[j]*blocksize+i; 92*360ee056SFande Kong } 93*360ee056SFande Kong ierr = MatSetValues(*interp,1,&row,nz,indices,values,INSERT_VALUES);CHKERRQ(ierr); 94*360ee056SFande Kong } 95*360ee056SFande Kong HYPRE_ParCSRMatrixRestoreRow(subinterp,subrow,(HYPRE_Int*)&nz,(HYPRE_Int**)&idx,&values); 96*360ee056SFande Kong } 97*360ee056SFande Kong ierr = PetscFree(indices);CHKERRQ(ierr); 98*360ee056SFande Kong ierr = MatAssemblyBegin(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 99*360ee056SFande Kong ierr = MatAssemblyEnd(*interp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 100*360ee056SFande Kong PetscFunctionReturn(0); 101*360ee056SFande Kong } 102*360ee056SFande Kong 103*360ee056SFande Kong PetscErrorCode PCSetUp_HMG(PC pc) 104*360ee056SFande Kong { 105*360ee056SFande Kong PetscErrorCode ierr; 106*360ee056SFande Kong Mat PA, submat; 107*360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 108*360ee056SFande Kong PC_HMG *hmg = (PC_HMG*) mg->innerctx; 109*360ee056SFande Kong MPI_Comm comm; 110*360ee056SFande Kong PetscInt level; 111*360ee056SFande Kong hypre_ParCSRMatrix **P_array, **A_array; 112*360ee056SFande Kong HYPRE_Solver hsolver; 113*360ee056SFande Kong hypre_ParAMGData *amg_data; 114*360ee056SFande Kong PetscInt num_levels; 115*360ee056SFande Kong PetscReal global_nonzeros, num_rows; 116*360ee056SFande Kong PetscReal *sparse; 117*360ee056SFande Kong 118*360ee056SFande Kong PetscFunctionBegin; 119*360ee056SFande Kong 120*360ee056SFande Kong ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 121*360ee056SFande Kong 122*360ee056SFande Kong if(pc->setupcalled){ 123*360ee056SFande Kong if (pc->flag == SAME_NONZERO_PATTERN && hmg->reuseinterp) { 124*360ee056SFande Kong ierr = PCMGSetGalerkin(pc,PC_MG_GALERKIN_PMAT);CHKERRQ(ierr); 125*360ee056SFande Kong ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 126*360ee056SFande Kong PetscFunctionReturn(0); 127*360ee056SFande Kong }else { 128*360ee056SFande Kong ierr = PCReset_MG(pc);CHKERRQ(ierr); 129*360ee056SFande Kong pc->setupcalled = PETSC_FALSE; 130*360ee056SFande Kong } 131*360ee056SFande Kong } 132*360ee056SFande Kong 133*360ee056SFande Kong if (!hmg->hypre){ 134*360ee056SFande Kong ierr = PCCreate(comm,&hmg->hypre);CHKERRQ(ierr); 135*360ee056SFande Kong ierr = PCSetType(hmg->hypre,PCHYPRE);CHKERRQ(ierr); 136*360ee056SFande Kong ierr = PCHYPRESetType(hmg->hypre,"boomeramg");CHKERRQ(ierr); 137*360ee056SFande Kong } 138*360ee056SFande Kong ierr = PCGetOperators(pc,NULL,&PA);CHKERRQ(ierr); 139*360ee056SFande Kong if(hmg->blocksize>1) { 140*360ee056SFande Kong ierr = PCHMGExtractSubMatrix_HMG(PA,&submat,MAT_INITIAL_MATRIX,hmg->blocksize);CHKERRQ(ierr); 141*360ee056SFande Kong PA = submat; 142*360ee056SFande Kong } 143*360ee056SFande Kong ierr = PCSetOperators(hmg->hypre,PA,PA);CHKERRQ(ierr); 144*360ee056SFande Kong if (hmg->blocksize>1){ 145*360ee056SFande Kong ierr = MatDestroy(&PA);CHKERRQ(ierr); 146*360ee056SFande Kong } 147*360ee056SFande Kong ierr = PCSetUseAmat(hmg->hypre,PETSC_FALSE);CHKERRQ(ierr); 148*360ee056SFande Kong ierr = PetscObjectOptionsBegin((PetscObject)hmg->hypre);CHKERRQ(ierr); 149*360ee056SFande Kong ierr = PCSetFromOptions_HYPRE(PetscOptionsObject,hmg->hypre);CHKERRQ(ierr); 150*360ee056SFande Kong ierr = PetscOptionsEnd();CHKERRQ(ierr); 151*360ee056SFande Kong ierr = PCSetUp(hmg->hypre); 152*360ee056SFande Kong 153*360ee056SFande Kong ierr = PCHYPREGetSolver(hmg->hypre,&hsolver);CHKERRQ(ierr); 154*360ee056SFande Kong amg_data = (hypre_ParAMGData*) (hsolver); 155*360ee056SFande Kong 156*360ee056SFande Kong if (!amg_data) SETERRQ(comm,PETSC_ERR_ARG_BADPTR,"Fails to setup Hypre \n"); 157*360ee056SFande Kong 158*360ee056SFande Kong num_levels = hypre_ParAMGDataNumLevels(amg_data); 159*360ee056SFande Kong ierr = PetscCalloc1(num_levels,&sparse);CHKERRQ(ierr); 160*360ee056SFande Kong 161*360ee056SFande Kong A_array= hypre_ParAMGDataAArray(amg_data); 162*360ee056SFande Kong P_array = hypre_ParAMGDataPArray(amg_data); 163*360ee056SFande Kong for(level=0;level<num_levels;level++){ 164*360ee056SFande Kong global_nonzeros = (PetscReal) hypre_ParCSRMatrixDNumNonzeros(A_array[level]); 165*360ee056SFande Kong num_rows = (PetscReal) hypre_ParCSRMatrixGlobalNumRows(A_array[level]); 166*360ee056SFande Kong sparse[num_levels-1-level] =global_nonzeros/(num_rows*num_rows); 167*360ee056SFande Kong } 168*360ee056SFande Kong ierr = PCMGSetLevels_MG(pc,num_levels,NULL);CHKERRQ(ierr); 169*360ee056SFande Kong ierr = PetscFree(sparse);CHKERRQ(ierr); 170*360ee056SFande Kong 171*360ee056SFande Kong for(level=num_levels-1;level>0;level--){ 172*360ee056SFande Kong Mat P=0, pmat=0; 173*360ee056SFande Kong Vec b, x,r; 174*360ee056SFande Kong if (hmg->blocksize>1){ 175*360ee056SFande Kong ierr = PCHMGExpandInterpolation_HMG(comm,P_array[num_levels-1-level],&P,hmg->blocksize);CHKERRQ(ierr); 176*360ee056SFande Kong }else{ 177*360ee056SFande Kong ierr = MatConvert_ParCSRMatrix_AIJ(comm, P_array[num_levels-1-level],MATAIJ, MAT_INPLACE_MATRIX,&P);CHKERRQ(ierr); 178*360ee056SFande Kong } 179*360ee056SFande Kong /*ierr = MatView(P,NULL);CHKERRQ(ierr);*/ 180*360ee056SFande Kong ierr = MatCreateVecs(P,&b,&r);CHKERRQ(ierr); 181*360ee056SFande Kong ierr = PCMGSetInterpolation(pc,level,P);CHKERRQ(ierr); 182*360ee056SFande Kong ierr = PCMGSetRestriction(pc,level,P);CHKERRQ(ierr); 183*360ee056SFande Kong ierr = MatDestroy(&P);CHKERRQ(ierr); 184*360ee056SFande Kong if ((level-1)>=0 && hmg->blocksize<=1) { 185*360ee056SFande Kong ierr = MatConvert_ParCSRMatrix_AIJ(comm, A_array[num_levels-level],MATAIJ, MAT_INPLACE_MATRIX,&pmat);CHKERRQ(ierr); 186*360ee056SFande Kong ierr = PCMGSetOperators(pc,level-1,pmat);CHKERRQ(ierr); 187*360ee056SFande Kong ierr = MatDestroy(&pmat);CHKERRQ(ierr); 188*360ee056SFande Kong } 189*360ee056SFande Kong ierr = PCMGSetRhs(pc,level-1,b);CHKERRQ(ierr); 190*360ee056SFande Kong 191*360ee056SFande Kong ierr = PCMGSetR(pc,level,r);CHKERRQ(ierr); 192*360ee056SFande Kong ierr = VecDestroy(&r);CHKERRQ(ierr); 193*360ee056SFande Kong 194*360ee056SFande Kong ierr = VecDuplicate(b,&x); 195*360ee056SFande Kong ierr = PCMGSetX(pc,level-1,x);CHKERRQ(ierr); 196*360ee056SFande Kong ierr = VecDestroy(&x);CHKERRQ(ierr); 197*360ee056SFande Kong ierr = VecDestroy(&b);CHKERRQ(ierr); 198*360ee056SFande Kong } 199*360ee056SFande Kong ierr = PCDestroy(&hmg->hypre);CHKERRQ(ierr); 200*360ee056SFande Kong hmg->hypre = 0; 201*360ee056SFande Kong ierr = PCMGSetGalerkin(pc,(hmg->blocksize>1) ? PC_MG_GALERKIN_PMAT:PC_MG_GALERKIN_NONE);CHKERRQ(ierr); 202*360ee056SFande Kong ierr = PCSetDM(pc,NULL);CHKERRQ(ierr); 203*360ee056SFande Kong ierr = PCSetUseAmat(pc,PETSC_FALSE);CHKERRQ(ierr); 204*360ee056SFande Kong ierr = PetscObjectOptionsBegin((PetscObject)pc);CHKERRQ(ierr); 205*360ee056SFande Kong ierr = PCSetFromOptions_MG(PetscOptionsObject,pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_HMG(), but cannot be called prior to PCMGSetLevels() */ 206*360ee056SFande Kong ierr = PetscOptionsEnd();CHKERRQ(ierr); 207*360ee056SFande Kong ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 208*360ee056SFande Kong PetscFunctionReturn(0); 209*360ee056SFande Kong } 210*360ee056SFande Kong 211*360ee056SFande Kong PetscErrorCode PCDestroy_HMG(PC pc) 212*360ee056SFande Kong { 213*360ee056SFande Kong PetscErrorCode ierr; 214*360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 215*360ee056SFande Kong PC_HMG *ctx = (PC_HMG*) mg->innerctx; 216*360ee056SFande Kong 217*360ee056SFande Kong PetscFunctionBegin; 218*360ee056SFande Kong ierr = PCDestroy(&ctx->hypre);CHKERRQ(ierr); 219*360ee056SFande Kong ierr = PetscFree(ctx);CHKERRQ(ierr); 220*360ee056SFande Kong ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 221*360ee056SFande Kong PetscFunctionReturn(0); 222*360ee056SFande Kong } 223*360ee056SFande Kong 224*360ee056SFande Kong PetscErrorCode PCView_HMG(PC pc,PetscViewer viewer) 225*360ee056SFande Kong { 226*360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 227*360ee056SFande Kong PC_HMG *ctx = (PC_HMG*) mg->innerctx; 228*360ee056SFande Kong PetscErrorCode ierr; 229*360ee056SFande Kong PetscBool iascii; 230*360ee056SFande Kong 231*360ee056SFande Kong PetscFunctionBegin; 232*360ee056SFande Kong ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 233*360ee056SFande Kong if (iascii) { 234*360ee056SFande Kong ierr = PetscViewerASCIIPrintf(viewer," Reuse interpolation %d\n",ctx->reuseinterp);CHKERRQ(ierr); 235*360ee056SFande Kong ierr = PetscViewerASCIIPrintf(viewer," Matrix block size %D \n",ctx->blocksize);CHKERRQ(ierr); 236*360ee056SFande Kong } 237*360ee056SFande Kong ierr = PCView_MG(pc,viewer);CHKERRQ(ierr); 238*360ee056SFande Kong PetscFunctionReturn(0); 239*360ee056SFande Kong } 240*360ee056SFande Kong 241*360ee056SFande Kong PetscErrorCode PCSetFromOptions_HMG(PetscOptionItems *PetscOptionsObject,PC pc) 242*360ee056SFande Kong { 243*360ee056SFande Kong PetscErrorCode ierr; 244*360ee056SFande Kong PC_MG *mg = (PC_MG*)pc->data; 245*360ee056SFande Kong PC_HMG *ctx = (PC_HMG*) mg->innerctx; 246*360ee056SFande Kong 247*360ee056SFande Kong PetscFunctionBegin; 248*360ee056SFande Kong ierr = PetscOptionsHead(PetscOptionsObject,"HMG");CHKERRQ(ierr); 249*360ee056SFande Kong ierr = PetscOptionsBool("-pc_hmg_reuse_interpolation","Reuse the interpolation operators when possible (cheaper, weaker when matrix entries change a lot)","None",ctx->reuseinterp,&ctx->reuseinterp,NULL);CHKERRQ(ierr); 250*360ee056SFande Kong ierr = PetscOptionsInt("-pc_hmg_pmat_blocksize","Block size for each grid point","hmg",ctx->blocksize,&ctx->blocksize,NULL);CHKERRQ(ierr); 251*360ee056SFande Kong ierr = PetscOptionsTail();CHKERRQ(ierr); 252*360ee056SFande Kong PetscFunctionReturn(0); 253*360ee056SFande Kong } 254*360ee056SFande Kong 255*360ee056SFande Kong /*MC 256*360ee056SFande Kong PCHMG - Hybrid of PETSc preconditioners (such as ASM, BJacobi, SOR, etc.) and Hypre BoomerAMG. BoomerAMG is used to coarsen matrix and generate 257*360ee056SFande Kong a sequence of coarse matrices and interpolations. The matrices and interpolations are employed to construct PCMG, and then any available 258*360ee056SFande Kong PETSc preconditioners can be chosen as smoothers and the coarse solver. 259*360ee056SFande Kong 260*360ee056SFande Kong Options Database Keys: 261*360ee056SFande Kong + -pc_hmg_reuse_interpolation <true | false> - Whether or not or not to reuse the interpolations. If true, it potentially save the compute time. 262*360ee056SFande Kong . -pc_hmg_pmat_blocksize - Block size of the underlying matrix. 263*360ee056SFande Kong 264*360ee056SFande Kong 265*360ee056SFande Kong Notes: 266*360ee056SFande Kong For multicomponent problems, we can just coarsen one submatrix associated with one particular component. In this way, the preconditioner setup 267*360ee056SFande Kong time is significantly reduced. One typical use case is neutron transport equations. There are many variables on each mesh vertex due to the 268*360ee056SFande Kong of angle and energy. Each variable, in fact, corresponds to the same PDEs but with different material properties. 269*360ee056SFande Kong 270*360ee056SFande Kong Level: beginner 271*360ee056SFande Kong 272*360ee056SFande Kong Concepts: additive Schwarz method 273*360ee056SFande Kong 274*360ee056SFande Kong References: 275*360ee056SFande 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 276*360ee056SFande Kong Newton-Krylov-Schwarz method with subspace-based coarsening and partition-based balancing for the multigroup neutron transport equations on 277*360ee056SFande Kong 3D unstructured meshes, arXiv preprint arXiv:1903.03659, 2019 278*360ee056SFande Kong 279*360ee056SFande Kong .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMG, PCHYPRE 280*360ee056SFande Kong 281*360ee056SFande Kong M*/ 282*360ee056SFande Kong PETSC_EXTERN PetscErrorCode PCCreate_HMG(PC pc) 283*360ee056SFande Kong { 284*360ee056SFande Kong PetscErrorCode ierr; 285*360ee056SFande Kong PC_HMG *hmg; 286*360ee056SFande Kong PC_MG *mg; 287*360ee056SFande Kong 288*360ee056SFande Kong PetscFunctionBegin; 289*360ee056SFande Kong /* if type was previously mg; must manually destroy it because call to PCSetType(pc,PCMG) will not destroy it */ 290*360ee056SFande Kong if (pc->ops->destroy) { 291*360ee056SFande Kong ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 292*360ee056SFande Kong pc->data = 0; 293*360ee056SFande Kong } 294*360ee056SFande Kong ierr = PetscFree(((PetscObject)pc)->type_name);CHKERRQ(ierr); 295*360ee056SFande Kong ((PetscObject)pc)->type_name = 0; 296*360ee056SFande Kong 297*360ee056SFande Kong ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); 298*360ee056SFande Kong ierr = PetscNew(&hmg);CHKERRQ(ierr); \ 299*360ee056SFande Kong 300*360ee056SFande Kong mg = (PC_MG*) pc->data; 301*360ee056SFande Kong mg->innerctx = hmg; 302*360ee056SFande Kong hmg->reuseinterp = PETSC_FALSE; 303*360ee056SFande Kong hmg->blocksize = 1; 304*360ee056SFande Kong 305*360ee056SFande Kong pc->ops->setfromoptions = PCSetFromOptions_HMG; 306*360ee056SFande Kong pc->ops->view = PCView_HMG; 307*360ee056SFande Kong pc->ops->destroy = PCDestroy_HMG; 308*360ee056SFande Kong pc->ops->setup = PCSetUp_HMG; 309*360ee056SFande Kong PetscFunctionReturn(0); 310*360ee056SFande Kong } 311