/* Creates hypre ijmatrix from PETSc matrix */ #include #include #include <../src/mat/impls/hypre/mhypre.h> #include <../src/mat/impls/aij/mpi/mpiaij.h> #include #include #include <_hypre_struct_mv.h> static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*); static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix); static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix); static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix); static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool); #undef __FUNCT__ #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate" static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij) { PetscErrorCode ierr; PetscInt i,n_d,n_o; const PetscInt *ia_d,*ia_o; PetscBool done_d=PETSC_FALSE,done_o=PETSC_FALSE; PetscInt *nnz_d=NULL,*nnz_o=NULL; PetscFunctionBegin; if (A_d) { /* determine number of nonzero entries in local diagonal part */ ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr); if (done_d) { ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr); for (i=0; irmap->rstart; rend = A->rmap->rend; cstart = A->cmap->rstart; cend = A->cmap->rend; PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij)); PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); { PetscBool same; Mat A_d,A_o; const PetscInt *colmap; ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr); if (same) { ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr); PetscFunctionReturn(0); } ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr); if (same) { ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr); PetscFunctionReturn(0); } ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr); if (same) { ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr); PetscFunctionReturn(0); } ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr); if (same) { ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr); PetscFunctionReturn(0); } } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatHYPRE_IJMatrixCopy" static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij) { PetscErrorCode ierr; PetscInt i,rstart,rend,ncols,nr,nc; const PetscScalar *values; const PetscInt *cols; PetscBool flg; PetscFunctionBegin; ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr); if (flg && nr == nc) { ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr); PetscFunctionReturn(0); } ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); if (flg) { ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr); PetscFunctionReturn(0); } PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); for (i=rstart; idata; hypre_ParCSRMatrix *par_matrix; hypre_AuxParCSRMatrix *aux_matrix; hypre_CSRMatrix *hdiag; PetscFunctionBegin; PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij); aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); hdiag = hypre_ParCSRMatrixDiag(par_matrix); /* this is the Hack part where we monkey directly with the hypre datastructures */ ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt)); ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt)); ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ" static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij) { PetscErrorCode ierr; Mat_MPIAIJ *pA = (Mat_MPIAIJ*)A->data; Mat_SeqAIJ *pdiag,*poffd; PetscInt i,*garray = pA->garray,*jj,cstart,*pjj; hypre_ParCSRMatrix *par_matrix; hypre_AuxParCSRMatrix *aux_matrix; hypre_CSRMatrix *hdiag,*hoffd; PetscFunctionBegin; pdiag = (Mat_SeqAIJ*) pA->A->data; poffd = (Mat_SeqAIJ*) pA->B->data; /* cstart is only valid for square MPIAIJ layed out in the usual way */ ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr); PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij); aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); hdiag = hypre_ParCSRMatrixDiag(par_matrix); hoffd = hypre_ParCSRMatrixOffd(par_matrix); /* this is the Hack part where we monkey directly with the hypre datastructures */ ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */ jj = (PetscInt*)hdiag->j; pjj = (PetscInt*)pdiag->j; for (i=0; inz; i++) jj[i] = cstart + pjj[i]; ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this If we hacked a hypre a bit more we might be able to avoid this step */ jj = (PetscInt*) hoffd->j; pjj = (PetscInt*) poffd->j; for (i=0; inz; i++) jj[i] = garray[pjj[i]]; ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar)); hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatConvert_AIJ_HYPRE" PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B) { Mat_HYPRE *hB; MPI_Comm comm = PetscObjectComm((PetscObject)A); PetscErrorCode ierr; PetscFunctionBegin; if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX"); if (reuse == MAT_REUSE_MATRIX) { /* always destroy the old matrix and create a new memory; hope this does not churn the memory too much. The problem is I do not know if it is possible to put the matrix back to its initial state so that we can directly copy the values the second time through. */ hB = (Mat_HYPRE*)((*B)->data); PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij)); } else { ierr = MatCreate(comm,B);CHKERRQ(ierr); ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr); ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); ierr = MatSetUp(*B);CHKERRQ(ierr); hB = (Mat_HYPRE*)((*B)->data); } ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr); ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr); ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatConvert_HYPRE_AIJ" PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B) { Mat_HYPRE *hA = (Mat_HYPRE*)A->data; hypre_ParCSRMatrix *parcsr; hypre_CSRMatrix *hdiag,*hoffd; MPI_Comm comm; PetscScalar *da,*oa,*aptr; PetscInt *dii,*djj,*oii,*ojj,*iptr; PetscInt i,dnnz,onnz,m,n; int type; PetscMPIInt size; PetscErrorCode ierr; PetscFunctionBegin; comm = PetscObjectComm((PetscObject)A); PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); if (reuse == MAT_REUSE_MATRIX) { PetscBool ismpiaij,isseqaij; ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported"); } ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); hdiag = hypre_ParCSRMatrixDiag(parcsr); hoffd = hypre_ParCSRMatrixOffd(parcsr); m = hypre_CSRMatrixNumRows(hdiag); n = hypre_CSRMatrixNumCols(hdiag); dnnz = hypre_CSRMatrixNumNonzeros(hdiag); onnz = hypre_CSRMatrixNumNonzeros(hoffd); if (reuse != MAT_REUSE_MATRIX) { ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr); ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr); ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr); } else { PetscInt nr; PetscBool done; if (size > 1) { Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows in diag part! %D != %D",nr,m); if (dii[nr] < dnnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros in diag part! reuse %D hypre %D",dii[nr],dnnz); ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr); } else { ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m); if (dii[nr] < dnnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros! reuse %D hypre %D",dii[nr],dnnz); ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr); } } ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr); ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr); iptr = djj; aptr = da; for (i=0; i 1) { PetscInt *offdj,*coffd; if (reuse != MAT_REUSE_MATRIX) { ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr); ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr); ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr); } else { Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); PetscInt nr,hr = hypre_CSRMatrixNumRows(hoffd); PetscBool done; ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr); if (nr != hr) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows in offdiag part! %D != %D",nr,hr); if (oii[nr] < onnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros in offdiag part! reuse %D hypre %D",oii[nr],onnz); ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr); } ierr = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); offdj = hypre_CSRMatrixJ(hoffd); coffd = hypre_ParCSRMatrixColMapOffd(parcsr); for (i=0; idata); d = (Mat_SeqAIJ*)b->A->data; o = (Mat_SeqAIJ*)b->B->data; d->free_a = PETSC_TRUE; d->free_ij = PETSC_TRUE; o->free_a = PETSC_TRUE; o->free_ij = PETSC_TRUE; } } else if (reuse != MAT_REUSE_MATRIX) { Mat_SeqAIJ* b; ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr); /* hack SeqAIJ */ b = (Mat_SeqAIJ*)((*B)->data); b->free_a = PETSC_TRUE; b->free_ij = PETSC_TRUE; } if (reuse == MAT_INPLACE_MATRIX) { ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); } PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMultTranspose_HYPRE" PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y) { PetscErrorCode ierr; PetscFunctionBegin; ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatMult_HYPRE" PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y) { PetscErrorCode ierr; PetscFunctionBegin; ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatHYPRE_MultKernel_Private" PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans) { Mat_HYPRE *hA = (Mat_HYPRE*)A->data; hypre_ParCSRMatrix *parcsr; hypre_ParVector *hx,*hy; PetscScalar *ax,*ay,*sax,*say; PetscErrorCode ierr; PetscFunctionBegin; PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx)); PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy)); ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); ierr = VecGetArray(y,&ay);CHKERRQ(ierr); if (trans) { HYPREReplacePointer(hA->x,ay,say); HYPREReplacePointer(hA->b,ax,sax); hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx); HYPREReplacePointer(hA->x,say,ay); HYPREReplacePointer(hA->b,sax,ax); } else { HYPREReplacePointer(hA->x,ax,sax); HYPREReplacePointer(hA->b,ay,say); hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy); HYPREReplacePointer(hA->x,sax,ax); HYPREReplacePointer(hA->b,say,ay); } ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatDestroy_HYPRE" PetscErrorCode MatDestroy_HYPRE(Mat A) { Mat_HYPRE *hA = (Mat_HYPRE*)A->data; PetscErrorCode ierr; PetscFunctionBegin; if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x)); if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b)); if (hA->ij) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij)); if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr);} ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr); ierr = PetscFree(A->data);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatSetUp_HYPRE" PetscErrorCode MatSetUp_HYPRE(Mat A) { Mat_HYPRE *hA = (Mat_HYPRE*)A->data; Vec x,b; PetscErrorCode ierr; PetscFunctionBegin; if (hA->x) PetscFunctionReturn(0); ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr); ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr); ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr); ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr); ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&b);CHKERRQ(ierr); PetscFunctionReturn(0); } #undef __FUNCT__ #define __FUNCT__ "MatAssemblyEnd_HYPRE" PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode) { Mat_HYPRE *hA = (Mat_HYPRE*)A->data; PetscErrorCode ierr; PetscFunctionBegin; PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij)); PetscFunctionReturn(0); } /* Is this needed? int type; hypre_ParCSRMatrix *parcsr; PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR supported"); PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); PetscStackCallStandard(hypre_MatvecCommPkgCreate,(parcsr)); */ #undef __FUNCT__ #define __FUNCT__ "MatCreate_HYPRE" PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B) { Mat_HYPRE *hB; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscNewLog(B,&hB);CHKERRQ(ierr); B->data = (void*)hB; B->rmap->bs = 1; B->assembled = PETSC_FALSE; B->ops->mult = MatMult_HYPRE; B->ops->multtranspose = MatMultTranspose_HYPRE; B->ops->setup = MatSetUp_HYPRE; B->ops->destroy = MatDestroy_HYPRE; B->ops->assemblyend = MatAssemblyEnd_HYPRE; ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr); ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr); ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr); PetscFunctionReturn(0); } #if 0 /* Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first which will corrupt the PETSc data structure if we did this. Need a work around to this problem. */ #include <_hypre_IJ_mv.h> #include #undef __FUNCT__ #define __FUNCT__ "MatHYPRE_IJMatrixLink" PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij) { PetscErrorCode ierr; PetscInt rstart,rend,cstart,cend; PetscBool flg; hypre_AuxParCSRMatrix *aux_matrix; PetscFunctionBegin; PetscValidHeaderSpecific(A,MAT_CLASSID,1); PetscValidType(A,1); PetscValidPointer(ij,2); ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices"); ierr = MatSetUp(A);CHKERRQ(ierr); ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr); rstart = A->rmap->rstart; rend = A->rmap->rend; cstart = A->cmap->rstart; cend = A->cmap->rend; PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij)); PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR)); PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij)); PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij)); hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; /* this is the Hack part where we monkey directly with the hypre datastructures */ PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij)); ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); } #endif