163c07aadSStefano Zampini 263c07aadSStefano Zampini /* 363c07aadSStefano Zampini Creates hypre ijmatrix from PETSc matrix 463c07aadSStefano Zampini */ 5225daaf8SStefano Zampini 6225daaf8SStefano Zampini /*MC 7225daaf8SStefano Zampini MATHYPRE - MATHYPRE = "hypre" - A matrix type to be used for sequential and parallel sparse matrices 8225daaf8SStefano Zampini based on the hypre IJ interface. 9225daaf8SStefano Zampini 10225daaf8SStefano Zampini Level: intermediate 11225daaf8SStefano Zampini 12225daaf8SStefano Zampini .seealso: MatCreate() 13225daaf8SStefano Zampini M*/ 14225daaf8SStefano Zampini 15dd9c0a25Sstefano_zampini #include <petscmathypre.h> 1663c07aadSStefano Zampini #include <petsc/private/matimpl.h> 1763c07aadSStefano Zampini #include <../src/mat/impls/hypre/mhypre.h> 1863c07aadSStefano Zampini #include <../src/mat/impls/aij/mpi/mpiaij.h> 1958968eb6SStefano Zampini #include <../src/vec/vec/impls/hypre/vhyp.h> 2058968eb6SStefano Zampini #include <HYPRE.h> 21c1a070e6SStefano Zampini #include <HYPRE_utilities.h> 22cd8bc7baSStefano Zampini #include <_hypre_parcsr_ls.h> 2363c07aadSStefano Zampini 2463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*); 2563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix); 2663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix); 2763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix); 2863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool); 29225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void*); 3063c07aadSStefano Zampini 3163c07aadSStefano Zampini #undef __FUNCT__ 3263c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate" 3363c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij) 3463c07aadSStefano Zampini { 3563c07aadSStefano Zampini PetscErrorCode ierr; 3663c07aadSStefano Zampini PetscInt i,n_d,n_o; 3763c07aadSStefano Zampini const PetscInt *ia_d,*ia_o; 3863c07aadSStefano Zampini PetscBool done_d=PETSC_FALSE,done_o=PETSC_FALSE; 3963c07aadSStefano Zampini PetscInt *nnz_d=NULL,*nnz_o=NULL; 4063c07aadSStefano Zampini 4163c07aadSStefano Zampini PetscFunctionBegin; 4263c07aadSStefano Zampini if (A_d) { /* determine number of nonzero entries in local diagonal part */ 4363c07aadSStefano Zampini ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr); 4463c07aadSStefano Zampini if (done_d) { 4563c07aadSStefano Zampini ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr); 4663c07aadSStefano Zampini for (i=0; i<n_d; i++) { 4763c07aadSStefano Zampini nnz_d[i] = ia_d[i+1] - ia_d[i]; 4863c07aadSStefano Zampini } 4963c07aadSStefano Zampini } 5063c07aadSStefano Zampini ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr); 5163c07aadSStefano Zampini } 5263c07aadSStefano Zampini if (A_o) { /* determine number of nonzero entries in local off-diagonal part */ 5363c07aadSStefano Zampini ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr); 5463c07aadSStefano Zampini if (done_o) { 5563c07aadSStefano Zampini ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr); 5663c07aadSStefano Zampini for (i=0; i<n_o; i++) { 5763c07aadSStefano Zampini nnz_o[i] = ia_o[i+1] - ia_o[i]; 5863c07aadSStefano Zampini } 5963c07aadSStefano Zampini } 6063c07aadSStefano Zampini ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr); 6163c07aadSStefano Zampini } 6263c07aadSStefano Zampini if (done_d) { /* set number of nonzeros in HYPRE IJ matrix */ 6363c07aadSStefano Zampini if (!done_o) { /* only diagonal part */ 6463c07aadSStefano Zampini ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr); 6563c07aadSStefano Zampini for (i=0; i<n_d; i++) { 6663c07aadSStefano Zampini nnz_o[i] = 0; 6763c07aadSStefano Zampini } 6863c07aadSStefano Zampini } 6963c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o)); 7063c07aadSStefano Zampini ierr = PetscFree(nnz_d);CHKERRQ(ierr); 7163c07aadSStefano Zampini ierr = PetscFree(nnz_o);CHKERRQ(ierr); 7263c07aadSStefano Zampini } 7363c07aadSStefano Zampini PetscFunctionReturn(0); 7463c07aadSStefano Zampini } 7563c07aadSStefano Zampini 7663c07aadSStefano Zampini #undef __FUNCT__ 7763c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_CreateFromMat" 7863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA) 7963c07aadSStefano Zampini { 8063c07aadSStefano Zampini PetscErrorCode ierr; 8163c07aadSStefano Zampini PetscInt rstart,rend,cstart,cend; 8263c07aadSStefano Zampini 8363c07aadSStefano Zampini PetscFunctionBegin; 84d975228cSstefano_zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 85d975228cSstefano_zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 8663c07aadSStefano Zampini rstart = A->rmap->rstart; 8763c07aadSStefano Zampini rend = A->rmap->rend; 8863c07aadSStefano Zampini cstart = A->cmap->rstart; 8963c07aadSStefano Zampini cend = A->cmap->rend; 9063c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij)); 9163c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); 9263c07aadSStefano Zampini { 9363c07aadSStefano Zampini PetscBool same; 9463c07aadSStefano Zampini Mat A_d,A_o; 9563c07aadSStefano Zampini const PetscInt *colmap; 9663c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr); 9763c07aadSStefano Zampini if (same) { 9863c07aadSStefano Zampini ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); 9963c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr); 10063c07aadSStefano Zampini PetscFunctionReturn(0); 10163c07aadSStefano Zampini } 10263c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr); 10363c07aadSStefano Zampini if (same) { 10463c07aadSStefano Zampini ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); 10563c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr); 10663c07aadSStefano Zampini PetscFunctionReturn(0); 10763c07aadSStefano Zampini } 10863c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr); 10963c07aadSStefano Zampini if (same) { 11063c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr); 11163c07aadSStefano Zampini PetscFunctionReturn(0); 11263c07aadSStefano Zampini } 11363c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr); 11463c07aadSStefano Zampini if (same) { 11563c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr); 11663c07aadSStefano Zampini PetscFunctionReturn(0); 11763c07aadSStefano Zampini } 11863c07aadSStefano Zampini } 11963c07aadSStefano Zampini PetscFunctionReturn(0); 12063c07aadSStefano Zampini } 12163c07aadSStefano Zampini 12263c07aadSStefano Zampini #undef __FUNCT__ 12363c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixCopy" 12463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij) 12563c07aadSStefano Zampini { 12663c07aadSStefano Zampini PetscErrorCode ierr; 12763c07aadSStefano Zampini PetscInt i,rstart,rend,ncols,nr,nc; 12863c07aadSStefano Zampini const PetscScalar *values; 12963c07aadSStefano Zampini const PetscInt *cols; 13063c07aadSStefano Zampini PetscBool flg; 13163c07aadSStefano Zampini 13263c07aadSStefano Zampini PetscFunctionBegin; 13363c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 13463c07aadSStefano Zampini ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr); 13563c07aadSStefano Zampini if (flg && nr == nc) { 13663c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr); 13763c07aadSStefano Zampini PetscFunctionReturn(0); 13863c07aadSStefano Zampini } 13963c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); 14063c07aadSStefano Zampini if (flg) { 14163c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr); 14263c07aadSStefano Zampini PetscFunctionReturn(0); 14363c07aadSStefano Zampini } 14463c07aadSStefano Zampini 14563c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 14663c07aadSStefano Zampini ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 14763c07aadSStefano Zampini for (i=rstart; i<rend; i++) { 14863c07aadSStefano Zampini ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 14963c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values)); 15063c07aadSStefano Zampini ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 15163c07aadSStefano Zampini } 15263c07aadSStefano Zampini PetscFunctionReturn(0); 15363c07aadSStefano Zampini } 15463c07aadSStefano Zampini 15563c07aadSStefano Zampini #undef __FUNCT__ 15663c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ" 15763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij) 15863c07aadSStefano Zampini { 15963c07aadSStefano Zampini PetscErrorCode ierr; 16063c07aadSStefano Zampini Mat_SeqAIJ *pdiag = (Mat_SeqAIJ*)A->data; 16158968eb6SStefano Zampini HYPRE_Int type; 16263c07aadSStefano Zampini hypre_ParCSRMatrix *par_matrix; 16363c07aadSStefano Zampini hypre_AuxParCSRMatrix *aux_matrix; 16463c07aadSStefano Zampini hypre_CSRMatrix *hdiag; 16563c07aadSStefano Zampini 16663c07aadSStefano Zampini PetscFunctionBegin; 16763c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 168ea9daf28SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type)); 169ea9daf28SStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 170ea9daf28SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix)); 17163c07aadSStefano Zampini hdiag = hypre_ParCSRMatrixDiag(par_matrix); 17263c07aadSStefano Zampini /* 17363c07aadSStefano Zampini this is the Hack part where we monkey directly with the hypre datastructures 17463c07aadSStefano Zampini */ 17563c07aadSStefano Zampini ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt)); 17663c07aadSStefano Zampini ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt)); 17763c07aadSStefano Zampini ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); 178ea9daf28SStefano Zampini 179ea9daf28SStefano Zampini aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); 18063c07aadSStefano Zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 18163c07aadSStefano Zampini PetscFunctionReturn(0); 18263c07aadSStefano Zampini } 18363c07aadSStefano Zampini 18463c07aadSStefano Zampini #undef __FUNCT__ 18563c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ" 18663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij) 18763c07aadSStefano Zampini { 18863c07aadSStefano Zampini PetscErrorCode ierr; 18963c07aadSStefano Zampini Mat_MPIAIJ *pA = (Mat_MPIAIJ*)A->data; 19063c07aadSStefano Zampini Mat_SeqAIJ *pdiag,*poffd; 19163c07aadSStefano Zampini PetscInt i,*garray = pA->garray,*jj,cstart,*pjj; 19258968eb6SStefano Zampini HYPRE_Int type; 19363c07aadSStefano Zampini hypre_ParCSRMatrix *par_matrix; 19463c07aadSStefano Zampini hypre_AuxParCSRMatrix *aux_matrix; 19563c07aadSStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 19663c07aadSStefano Zampini 19763c07aadSStefano Zampini PetscFunctionBegin; 19863c07aadSStefano Zampini pdiag = (Mat_SeqAIJ*) pA->A->data; 19963c07aadSStefano Zampini poffd = (Mat_SeqAIJ*) pA->B->data; 20063c07aadSStefano Zampini /* cstart is only valid for square MPIAIJ layed out in the usual way */ 20163c07aadSStefano Zampini ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr); 20263c07aadSStefano Zampini 20363c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 204ea9daf28SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type)); 205ea9daf28SStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 206ea9daf28SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix)); 20763c07aadSStefano Zampini hdiag = hypre_ParCSRMatrixDiag(par_matrix); 20863c07aadSStefano Zampini hoffd = hypre_ParCSRMatrixOffd(par_matrix); 20963c07aadSStefano Zampini 21063c07aadSStefano Zampini /* 21163c07aadSStefano Zampini this is the Hack part where we monkey directly with the hypre datastructures 21263c07aadSStefano Zampini */ 21363c07aadSStefano Zampini ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); 21463c07aadSStefano Zampini /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */ 21563c07aadSStefano Zampini jj = (PetscInt*)hdiag->j; 21663c07aadSStefano Zampini pjj = (PetscInt*)pdiag->j; 21763c07aadSStefano Zampini for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i]; 21863c07aadSStefano Zampini ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); 21963c07aadSStefano Zampini ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); 22063c07aadSStefano Zampini /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this 22163c07aadSStefano Zampini If we hacked a hypre a bit more we might be able to avoid this step */ 22263c07aadSStefano Zampini jj = (PetscInt*) hoffd->j; 22363c07aadSStefano Zampini pjj = (PetscInt*) poffd->j; 22463c07aadSStefano Zampini for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]]; 22563c07aadSStefano Zampini ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar)); 22663c07aadSStefano Zampini 227ea9daf28SStefano Zampini aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); 22863c07aadSStefano Zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 22963c07aadSStefano Zampini PetscFunctionReturn(0); 23063c07aadSStefano Zampini } 23163c07aadSStefano Zampini 23263c07aadSStefano Zampini #undef __FUNCT__ 2332df22349SStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_IS" 2342df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B) 2352df22349SStefano Zampini { 2362df22349SStefano Zampini Mat_HYPRE* mhA = (Mat_HYPRE*)(A->data); 2372df22349SStefano Zampini Mat lA; 2382df22349SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 2392df22349SStefano Zampini IS is; 2402df22349SStefano Zampini hypre_ParCSRMatrix *hA; 2412df22349SStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 2422df22349SStefano Zampini MPI_Comm comm; 2432df22349SStefano Zampini PetscScalar *hdd,*hod,*aa,*data; 2447d968826Sstefano_zampini HYPRE_Int *col_map_offd,*hdi,*hdj,*hoi,*hoj; 2452df22349SStefano Zampini PetscInt *ii,*jj,*iptr,*jptr; 2462df22349SStefano Zampini PetscInt cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N; 24758968eb6SStefano Zampini HYPRE_Int type; 2482df22349SStefano Zampini PetscErrorCode ierr; 2492df22349SStefano Zampini 2502df22349SStefano Zampini PetscFunctionBegin; 251a1787963SStefano Zampini comm = PetscObjectComm((PetscObject)A); 2522df22349SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type)); 2532df22349SStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 2542df22349SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA)); 2552df22349SStefano Zampini M = hypre_ParCSRMatrixGlobalNumRows(hA); 2562df22349SStefano Zampini N = hypre_ParCSRMatrixGlobalNumCols(hA); 2572df22349SStefano Zampini str = hypre_ParCSRMatrixFirstRowIndex(hA); 2582df22349SStefano Zampini stc = hypre_ParCSRMatrixFirstColDiag(hA); 2592df22349SStefano Zampini hdiag = hypre_ParCSRMatrixDiag(hA); 2602df22349SStefano Zampini hoffd = hypre_ParCSRMatrixOffd(hA); 2612df22349SStefano Zampini dr = hypre_CSRMatrixNumRows(hdiag); 2622df22349SStefano Zampini dc = hypre_CSRMatrixNumCols(hdiag); 2632df22349SStefano Zampini nnz = hypre_CSRMatrixNumNonzeros(hdiag); 2642df22349SStefano Zampini hdi = hypre_CSRMatrixI(hdiag); 2652df22349SStefano Zampini hdj = hypre_CSRMatrixJ(hdiag); 2662df22349SStefano Zampini hdd = hypre_CSRMatrixData(hdiag); 2672df22349SStefano Zampini oc = hypre_CSRMatrixNumCols(hoffd); 2682df22349SStefano Zampini nnz += hypre_CSRMatrixNumNonzeros(hoffd); 2692df22349SStefano Zampini hoi = hypre_CSRMatrixI(hoffd); 2702df22349SStefano Zampini hoj = hypre_CSRMatrixJ(hoffd); 2712df22349SStefano Zampini hod = hypre_CSRMatrixData(hoffd); 2722df22349SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 2732df22349SStefano Zampini PetscInt *aux; 2742df22349SStefano Zampini 2752df22349SStefano Zampini /* generate l2g maps for rows and cols */ 2762df22349SStefano Zampini ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr); 2772df22349SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 2782df22349SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 2792df22349SStefano Zampini col_map_offd = hypre_ParCSRMatrixColMapOffd(hA); 2802df22349SStefano Zampini ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr); 2812df22349SStefano Zampini for (i=0; i<dc; i++) aux[i] = i+stc; 2822df22349SStefano Zampini for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i]; 2832df22349SStefano Zampini ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 2842df22349SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 2852df22349SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 2862df22349SStefano Zampini /* create MATIS object */ 2872df22349SStefano Zampini ierr = MatCreate(comm,B);CHKERRQ(ierr); 2882df22349SStefano Zampini ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr); 2892df22349SStefano Zampini ierr = MatSetType(*B,MATIS);CHKERRQ(ierr); 2902df22349SStefano Zampini ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr); 2912df22349SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 2922df22349SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 2932df22349SStefano Zampini 2942df22349SStefano Zampini /* allocate CSR for local matrix */ 2952df22349SStefano Zampini ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr); 2962df22349SStefano Zampini ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr); 2972df22349SStefano Zampini ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr); 2982df22349SStefano Zampini } else { 2992df22349SStefano Zampini PetscInt nr; 3002df22349SStefano Zampini PetscBool done; 3012df22349SStefano Zampini ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr); 3022df22349SStefano Zampini ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr); 3032df22349SStefano Zampini if (nr != dr) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of rows in local mat! %D != %D",nr,dr); 3042df22349SStefano Zampini if (iptr[nr] < nnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros in local mat! reuse %D requested %D",iptr[nr],nnz); 3052df22349SStefano Zampini ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr); 3062df22349SStefano Zampini } 3072df22349SStefano Zampini /* merge local matrices */ 3082df22349SStefano Zampini ii = iptr; 3092df22349SStefano Zampini jj = jptr; 3102df22349SStefano Zampini aa = data; 3112df22349SStefano Zampini *ii = *(hdi++) + *(hoi++); 3122df22349SStefano Zampini for (jd=0,jo=0,cum=0; *ii<nnz; cum++) { 3132df22349SStefano Zampini PetscScalar *aold = aa; 3142df22349SStefano Zampini PetscInt *jold = jj,nc = jd+jo; 3152df22349SStefano Zampini for (; jd<*hdi; jd++) { *jj++ = *hdj++; *aa++ = *hdd++; } 3162df22349SStefano Zampini for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; } 3172df22349SStefano Zampini *(++ii) = *(hdi++) + *(hoi++); 3182df22349SStefano Zampini ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr); 3192df22349SStefano Zampini } 3202df22349SStefano Zampini for (; cum<dr; cum++) *(++ii) = nnz; 3212df22349SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 322a033916dSStefano Zampini Mat_SeqAIJ* a; 323a033916dSStefano Zampini 3242df22349SStefano Zampini ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr); 3252df22349SStefano Zampini ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr); 326a033916dSStefano Zampini /* hack SeqAIJ */ 327a033916dSStefano Zampini a = (Mat_SeqAIJ*)(lA->data); 328a033916dSStefano Zampini a->free_a = PETSC_TRUE; 329a033916dSStefano Zampini a->free_ij = PETSC_TRUE; 3302df22349SStefano Zampini ierr = MatDestroy(&lA);CHKERRQ(ierr); 3312df22349SStefano Zampini } 3322df22349SStefano Zampini ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3332df22349SStefano Zampini ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3342df22349SStefano Zampini if (reuse == MAT_INPLACE_MATRIX) { 3352df22349SStefano Zampini ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 3362df22349SStefano Zampini } 3372df22349SStefano Zampini PetscFunctionReturn(0); 3382df22349SStefano Zampini } 3392df22349SStefano Zampini 3402df22349SStefano Zampini #undef __FUNCT__ 34163c07aadSStefano Zampini #define __FUNCT__ "MatConvert_AIJ_HYPRE" 34263c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B) 34363c07aadSStefano Zampini { 34463c07aadSStefano Zampini Mat_HYPRE *hB; 34563c07aadSStefano Zampini MPI_Comm comm = PetscObjectComm((PetscObject)A); 34663c07aadSStefano Zampini PetscErrorCode ierr; 34763c07aadSStefano Zampini 34863c07aadSStefano Zampini PetscFunctionBegin; 34963c07aadSStefano Zampini if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX"); 35063c07aadSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 35163c07aadSStefano Zampini /* always destroy the old matrix and create a new memory; 35263c07aadSStefano Zampini hope this does not churn the memory too much. The problem 35363c07aadSStefano Zampini is I do not know if it is possible to put the matrix back to 35463c07aadSStefano Zampini its initial state so that we can directly copy the values 35563c07aadSStefano Zampini the second time through. */ 35663c07aadSStefano Zampini hB = (Mat_HYPRE*)((*B)->data); 35763c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij)); 35863c07aadSStefano Zampini } else { 35963c07aadSStefano Zampini ierr = MatCreate(comm,B);CHKERRQ(ierr); 36063c07aadSStefano Zampini ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr); 36163c07aadSStefano Zampini ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 36263c07aadSStefano Zampini hB = (Mat_HYPRE*)((*B)->data); 36363c07aadSStefano Zampini } 36463c07aadSStefano Zampini ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr); 36563c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr); 3664ec6421dSstefano_zampini (*B)->preallocated = PETSC_TRUE; 36763c07aadSStefano Zampini ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 36863c07aadSStefano Zampini ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 36963c07aadSStefano Zampini PetscFunctionReturn(0); 37063c07aadSStefano Zampini } 37163c07aadSStefano Zampini 37263c07aadSStefano Zampini #undef __FUNCT__ 37363c07aadSStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_AIJ" 374ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B) 37563c07aadSStefano Zampini { 37663c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 37763c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 37863c07aadSStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 37963c07aadSStefano Zampini MPI_Comm comm; 38063c07aadSStefano Zampini PetscScalar *da,*oa,*aptr; 38163c07aadSStefano Zampini PetscInt *dii,*djj,*oii,*ojj,*iptr; 38263c07aadSStefano Zampini PetscInt i,dnnz,onnz,m,n; 38358968eb6SStefano Zampini HYPRE_Int type; 38463c07aadSStefano Zampini PetscMPIInt size; 38563c07aadSStefano Zampini PetscErrorCode ierr; 38663c07aadSStefano Zampini 38763c07aadSStefano Zampini PetscFunctionBegin; 38863c07aadSStefano Zampini comm = PetscObjectComm((PetscObject)A); 38963c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 39063c07aadSStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 39163c07aadSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 39263c07aadSStefano Zampini PetscBool ismpiaij,isseqaij; 39363c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 39463c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 39563c07aadSStefano Zampini if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported"); 39663c07aadSStefano Zampini } 39763c07aadSStefano Zampini ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 39863c07aadSStefano Zampini 39963c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 40063c07aadSStefano Zampini hdiag = hypre_ParCSRMatrixDiag(parcsr); 40163c07aadSStefano Zampini hoffd = hypre_ParCSRMatrixOffd(parcsr); 40263c07aadSStefano Zampini m = hypre_CSRMatrixNumRows(hdiag); 40363c07aadSStefano Zampini n = hypre_CSRMatrixNumCols(hdiag); 40463c07aadSStefano Zampini dnnz = hypre_CSRMatrixNumNonzeros(hdiag); 40563c07aadSStefano Zampini onnz = hypre_CSRMatrixNumNonzeros(hoffd); 406225daaf8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 40763c07aadSStefano Zampini ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr); 40863c07aadSStefano Zampini ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr); 40963c07aadSStefano Zampini ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr); 410225daaf8SStefano Zampini } else if (reuse == MAT_REUSE_MATRIX) { 41163c07aadSStefano Zampini PetscInt nr; 41263c07aadSStefano Zampini PetscBool done; 41363c07aadSStefano Zampini if (size > 1) { 41463c07aadSStefano Zampini Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); 41563c07aadSStefano Zampini 41663c07aadSStefano Zampini ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); 41763c07aadSStefano Zampini 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); 41863c07aadSStefano Zampini 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); 41963c07aadSStefano Zampini ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr); 42063c07aadSStefano Zampini } else { 42163c07aadSStefano Zampini ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); 42263c07aadSStefano Zampini if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m); 42363c07aadSStefano Zampini 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); 42463c07aadSStefano Zampini ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr); 42563c07aadSStefano Zampini } 426225daaf8SStefano Zampini } else { /* MAT_INPLACE_MATRIX */ 4277d968826Sstefano_zampini dii = (PetscInt*)hypre_CSRMatrixI(hdiag); 4287d968826Sstefano_zampini djj = (PetscInt*)hypre_CSRMatrixJ(hdiag); 429225daaf8SStefano Zampini da = hypre_CSRMatrixData(hdiag); 43063c07aadSStefano Zampini } 43163c07aadSStefano Zampini ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 43263c07aadSStefano Zampini ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr); 43363c07aadSStefano Zampini ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr); 43463c07aadSStefano Zampini iptr = djj; 43563c07aadSStefano Zampini aptr = da; 43663c07aadSStefano Zampini for (i=0; i<m; i++) { 43763c07aadSStefano Zampini PetscInt nc = dii[i+1]-dii[i]; 43863c07aadSStefano Zampini ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr); 43963c07aadSStefano Zampini iptr += nc; 44063c07aadSStefano Zampini aptr += nc; 44163c07aadSStefano Zampini } 44263c07aadSStefano Zampini if (size > 1) { 4437d968826Sstefano_zampini HYPRE_Int *offdj,*coffd; 44463c07aadSStefano Zampini 445225daaf8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 44663c07aadSStefano Zampini ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr); 44763c07aadSStefano Zampini ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr); 44863c07aadSStefano Zampini ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr); 449225daaf8SStefano Zampini } else if (reuse == MAT_REUSE_MATRIX) { 45063c07aadSStefano Zampini Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); 45163c07aadSStefano Zampini PetscInt nr,hr = hypre_CSRMatrixNumRows(hoffd); 45263c07aadSStefano Zampini PetscBool done; 45363c07aadSStefano Zampini 45463c07aadSStefano Zampini ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr); 45563c07aadSStefano Zampini 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); 45663c07aadSStefano Zampini 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); 45763c07aadSStefano Zampini ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr); 458225daaf8SStefano Zampini } else { /* MAT_INPLACE_MATRIX */ 4597d968826Sstefano_zampini oii = (PetscInt*)hypre_CSRMatrixI(hoffd); 4607d968826Sstefano_zampini ojj = (PetscInt*)hypre_CSRMatrixJ(hoffd); 461225daaf8SStefano Zampini oa = hypre_CSRMatrixData(hoffd); 46263c07aadSStefano Zampini } 46363c07aadSStefano Zampini ierr = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 46463c07aadSStefano Zampini offdj = hypre_CSRMatrixJ(hoffd); 46563c07aadSStefano Zampini coffd = hypre_ParCSRMatrixColMapOffd(parcsr); 46663c07aadSStefano Zampini for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]]; 46763c07aadSStefano Zampini ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr); 46863c07aadSStefano Zampini iptr = ojj; 46963c07aadSStefano Zampini aptr = oa; 47063c07aadSStefano Zampini for (i=0; i<m; i++) { 47163c07aadSStefano Zampini PetscInt nc = oii[i+1]-oii[i]; 47263c07aadSStefano Zampini ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr); 47363c07aadSStefano Zampini iptr += nc; 47463c07aadSStefano Zampini aptr += nc; 47563c07aadSStefano Zampini } 476225daaf8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 47763c07aadSStefano Zampini Mat_MPIAIJ *b; 47863c07aadSStefano Zampini Mat_SeqAIJ *d,*o; 479225daaf8SStefano Zampini 48041073d71Sstefano_zampini ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,B);CHKERRQ(ierr); 48163c07aadSStefano Zampini /* hack MPIAIJ */ 48263c07aadSStefano Zampini b = (Mat_MPIAIJ*)((*B)->data); 48363c07aadSStefano Zampini d = (Mat_SeqAIJ*)b->A->data; 48463c07aadSStefano Zampini o = (Mat_SeqAIJ*)b->B->data; 48563c07aadSStefano Zampini d->free_a = PETSC_TRUE; 48663c07aadSStefano Zampini d->free_ij = PETSC_TRUE; 48763c07aadSStefano Zampini o->free_a = PETSC_TRUE; 48863c07aadSStefano Zampini o->free_ij = PETSC_TRUE; 489225daaf8SStefano Zampini } else if (reuse == MAT_INPLACE_MATRIX) { 490225daaf8SStefano Zampini Mat T; 49141073d71Sstefano_zampini ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,&T);CHKERRQ(ierr); 492225daaf8SStefano Zampini hypre_CSRMatrixI(hdiag) = NULL; 493225daaf8SStefano Zampini hypre_CSRMatrixJ(hdiag) = NULL; 494225daaf8SStefano Zampini hypre_CSRMatrixData(hdiag) = NULL; 495225daaf8SStefano Zampini hypre_CSRMatrixI(hoffd) = NULL; 496225daaf8SStefano Zampini hypre_CSRMatrixJ(hoffd) = NULL; 497225daaf8SStefano Zampini hypre_CSRMatrixData(hoffd) = NULL; 498225daaf8SStefano Zampini ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr); 49963c07aadSStefano Zampini } 500225daaf8SStefano Zampini } else { 501225daaf8SStefano Zampini oii = NULL; 502225daaf8SStefano Zampini ojj = NULL; 503225daaf8SStefano Zampini oa = NULL; 504225daaf8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 50563c07aadSStefano Zampini Mat_SeqAIJ* b; 50663c07aadSStefano Zampini ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr); 50763c07aadSStefano Zampini /* hack SeqAIJ */ 50863c07aadSStefano Zampini b = (Mat_SeqAIJ*)((*B)->data); 50963c07aadSStefano Zampini b->free_a = PETSC_TRUE; 51063c07aadSStefano Zampini b->free_ij = PETSC_TRUE; 511225daaf8SStefano Zampini } else if (reuse == MAT_INPLACE_MATRIX) { 512225daaf8SStefano Zampini Mat T; 513225daaf8SStefano Zampini ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,&T);CHKERRQ(ierr); 514225daaf8SStefano Zampini hypre_CSRMatrixI(hdiag) = NULL; 515225daaf8SStefano Zampini hypre_CSRMatrixJ(hdiag) = NULL; 516225daaf8SStefano Zampini hypre_CSRMatrixData(hdiag) = NULL; 517225daaf8SStefano Zampini ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr); 51863c07aadSStefano Zampini } 519225daaf8SStefano Zampini } 520225daaf8SStefano Zampini 521225daaf8SStefano Zampini /* we have to use hypre_Tfree to free the arrays */ 52263c07aadSStefano Zampini if (reuse == MAT_INPLACE_MATRIX) { 523225daaf8SStefano Zampini void *ptrs[6] = {dii,djj,da,oii,ojj,oa}; 524225daaf8SStefano Zampini const char *names[6] = {"_hypre_csr_dii", 525225daaf8SStefano Zampini "_hypre_csr_djj", 526225daaf8SStefano Zampini "_hypre_csr_da", 527225daaf8SStefano Zampini "_hypre_csr_oii", 528225daaf8SStefano Zampini "_hypre_csr_ojj", 529225daaf8SStefano Zampini "_hypre_csr_oa"}; 530225daaf8SStefano Zampini for (i=0; i<6; i++) { 531225daaf8SStefano Zampini PetscContainer c; 532225daaf8SStefano Zampini 533225daaf8SStefano Zampini ierr = PetscContainerCreate(comm,&c);CHKERRQ(ierr); 534225daaf8SStefano Zampini ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr); 535225daaf8SStefano Zampini ierr = PetscContainerSetUserDestroy(c,hypre_array_destroy);CHKERRQ(ierr); 536225daaf8SStefano Zampini ierr = PetscObjectCompose((PetscObject)(*B),names[i],(PetscObject)c);CHKERRQ(ierr); 537225daaf8SStefano Zampini ierr = PetscContainerDestroy(&c);CHKERRQ(ierr); 538225daaf8SStefano Zampini } 53963c07aadSStefano Zampini } 54063c07aadSStefano Zampini PetscFunctionReturn(0); 54163c07aadSStefano Zampini } 54263c07aadSStefano Zampini 54363c07aadSStefano Zampini #undef __FUNCT__ 544613e5ff0Sstefano_zampini #define __FUNCT__ "MatAIJGetParCSR_Private" 545613e5ff0Sstefano_zampini static PetscErrorCode MatAIJGetParCSR_Private(Mat A, hypre_ParCSRMatrix **hA) 546c1a070e6SStefano Zampini { 547613e5ff0Sstefano_zampini hypre_ParCSRMatrix *tA; 548c1a070e6SStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 549c1a070e6SStefano Zampini Mat_SeqAIJ *diag,*offd; 550c1a070e6SStefano Zampini PetscInt *garray,noffd,dnnz,onnz,*row_starts,*col_starts; 551c1a070e6SStefano Zampini MPI_Comm comm = PetscObjectComm((PetscObject)A); 552613e5ff0Sstefano_zampini PetscBool ismpiaij,isseqaij; 553c1a070e6SStefano Zampini PetscErrorCode ierr; 554c1a070e6SStefano Zampini 555c1a070e6SStefano Zampini PetscFunctionBegin; 556c1a070e6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 557c1a070e6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 558c1a070e6SStefano Zampini if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type); 559c1a070e6SStefano Zampini if (ismpiaij) { 560c1a070e6SStefano Zampini Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data); 561c1a070e6SStefano Zampini 562c1a070e6SStefano Zampini diag = (Mat_SeqAIJ*)a->A->data; 563c1a070e6SStefano Zampini offd = (Mat_SeqAIJ*)a->B->data; 564c1a070e6SStefano Zampini garray = a->garray; 565c1a070e6SStefano Zampini noffd = a->B->cmap->N; 566c1a070e6SStefano Zampini dnnz = diag->nz; 567c1a070e6SStefano Zampini onnz = offd->nz; 568c1a070e6SStefano Zampini } else { 569c1a070e6SStefano Zampini diag = (Mat_SeqAIJ*)A->data; 570c1a070e6SStefano Zampini offd = NULL; 571c1a070e6SStefano Zampini garray = NULL; 572c1a070e6SStefano Zampini noffd = 0; 573c1a070e6SStefano Zampini dnnz = diag->nz; 574c1a070e6SStefano Zampini onnz = 0; 575c1a070e6SStefano Zampini } 576225daaf8SStefano Zampini 577c1a070e6SStefano Zampini /* create a temporary ParCSR */ 578c1a070e6SStefano Zampini if (HYPRE_AssumedPartitionCheck()) { 579c1a070e6SStefano Zampini PetscMPIInt myid; 580c1a070e6SStefano Zampini 581c1a070e6SStefano Zampini ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr); 582c1a070e6SStefano Zampini row_starts = A->rmap->range + myid; 583c1a070e6SStefano Zampini col_starts = A->cmap->range + myid; 584c1a070e6SStefano Zampini } else { 585c1a070e6SStefano Zampini row_starts = A->rmap->range; 586c1a070e6SStefano Zampini col_starts = A->cmap->range; 587c1a070e6SStefano Zampini } 5887d968826Sstefano_zampini tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,(HYPRE_Int*)row_starts,(HYPRE_Int*)col_starts,noffd,dnnz,onnz); 589c1a070e6SStefano Zampini hypre_ParCSRMatrixSetRowStartsOwner(tA,0); 590c1a070e6SStefano Zampini hypre_ParCSRMatrixSetColStartsOwner(tA,0); 591c1a070e6SStefano Zampini 592225daaf8SStefano Zampini /* set diagonal part */ 593c1a070e6SStefano Zampini hdiag = hypre_ParCSRMatrixDiag(tA); 5947d968826Sstefano_zampini hypre_CSRMatrixI(hdiag) = (HYPRE_Int*)diag->i; 5957d968826Sstefano_zampini hypre_CSRMatrixJ(hdiag) = (HYPRE_Int*)diag->j; 596c1a070e6SStefano Zampini hypre_CSRMatrixData(hdiag) = diag->a; 597c1a070e6SStefano Zampini hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz; 598c1a070e6SStefano Zampini hypre_CSRMatrixSetRownnz(hdiag); 599c1a070e6SStefano Zampini hypre_CSRMatrixSetDataOwner(hdiag,0); 600c1a070e6SStefano Zampini 601225daaf8SStefano Zampini /* set offdiagonal part */ 602c1a070e6SStefano Zampini hoffd = hypre_ParCSRMatrixOffd(tA); 603c1a070e6SStefano Zampini if (offd) { 6047d968826Sstefano_zampini hypre_CSRMatrixI(hoffd) = (HYPRE_Int*)offd->i; 6057d968826Sstefano_zampini hypre_CSRMatrixJ(hoffd) = (HYPRE_Int*)offd->j; 606c1a070e6SStefano Zampini hypre_CSRMatrixData(hoffd) = offd->a; 607c1a070e6SStefano Zampini hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz; 608c1a070e6SStefano Zampini hypre_CSRMatrixSetRownnz(hoffd); 609c1a070e6SStefano Zampini hypre_CSRMatrixSetDataOwner(hoffd,0); 610c1a070e6SStefano Zampini hypre_ParCSRMatrixSetNumNonzeros(tA); 6117d968826Sstefano_zampini hypre_ParCSRMatrixColMapOffd(tA) = (HYPRE_Int*)garray; 612c1a070e6SStefano Zampini } 613613e5ff0Sstefano_zampini *hA = tA; 614613e5ff0Sstefano_zampini PetscFunctionReturn(0); 615613e5ff0Sstefano_zampini } 616c1a070e6SStefano Zampini 617613e5ff0Sstefano_zampini #undef __FUNCT__ 618613e5ff0Sstefano_zampini #define __FUNCT__ "MatAIJRestoreParCSR_Private" 619613e5ff0Sstefano_zampini static PetscErrorCode MatAIJRestoreParCSR_Private(Mat A, hypre_ParCSRMatrix **hA) 620613e5ff0Sstefano_zampini { 621613e5ff0Sstefano_zampini hypre_CSRMatrix *hdiag,*hoffd; 622c1a070e6SStefano Zampini 623613e5ff0Sstefano_zampini PetscFunctionBegin; 624613e5ff0Sstefano_zampini hdiag = hypre_ParCSRMatrixDiag(*hA); 625613e5ff0Sstefano_zampini hoffd = hypre_ParCSRMatrixOffd(*hA); 626c1a070e6SStefano Zampini /* set pointers to NULL before destroying tA */ 627c1a070e6SStefano Zampini hypre_CSRMatrixI(hdiag) = NULL; 628c1a070e6SStefano Zampini hypre_CSRMatrixJ(hdiag) = NULL; 629c1a070e6SStefano Zampini hypre_CSRMatrixData(hdiag) = NULL; 630c1a070e6SStefano Zampini hypre_CSRMatrixI(hoffd) = NULL; 631c1a070e6SStefano Zampini hypre_CSRMatrixJ(hoffd) = NULL; 632c1a070e6SStefano Zampini hypre_CSRMatrixData(hoffd) = NULL; 633613e5ff0Sstefano_zampini hypre_ParCSRMatrixColMapOffd(*hA) = NULL; 634613e5ff0Sstefano_zampini hypre_ParCSRMatrixDestroy(*hA); 635613e5ff0Sstefano_zampini *hA = NULL; 636613e5ff0Sstefano_zampini PetscFunctionReturn(0); 637613e5ff0Sstefano_zampini } 638613e5ff0Sstefano_zampini 639613e5ff0Sstefano_zampini /* calls RAP from BoomerAMG: 640613e5ff0Sstefano_zampini the resulting ParCSR will not own the column and row starts */ 641613e5ff0Sstefano_zampini #undef __FUNCT__ 642613e5ff0Sstefano_zampini #define __FUNCT__ "MatHYPRE_ParCSR_RAP" 643613e5ff0Sstefano_zampini static PetscErrorCode MatHYPRE_ParCSR_RAP(hypre_ParCSRMatrix *hR, hypre_ParCSRMatrix *hA, 644613e5ff0Sstefano_zampini hypre_ParCSRMatrix *hP, hypre_ParCSRMatrix **hRAP) 645613e5ff0Sstefano_zampini { 646613e5ff0Sstefano_zampini PetscErrorCode ierr; 647613e5ff0Sstefano_zampini HYPRE_Int P_owns_col_starts,R_owns_row_starts; 648613e5ff0Sstefano_zampini 649613e5ff0Sstefano_zampini PetscFunctionBegin; 650613e5ff0Sstefano_zampini P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(hP); 651613e5ff0Sstefano_zampini R_owns_row_starts = hypre_ParCSRMatrixOwnsRowStarts(hR); 652613e5ff0Sstefano_zampini PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(hR,hA,hP,hRAP)); 653613e5ff0Sstefano_zampini PetscStackCallStandard(hypre_ParCSRMatrixSetNumNonzeros,(*hRAP)); 654613e5ff0Sstefano_zampini /* hypre_BoomerAMGBuildCoarseOperator steals the col_starts from P and the row_starts from R */ 655613e5ff0Sstefano_zampini hypre_ParCSRMatrixSetRowStartsOwner(*hRAP,0); 656613e5ff0Sstefano_zampini hypre_ParCSRMatrixSetColStartsOwner(*hRAP,0); 657613e5ff0Sstefano_zampini if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(hP,1); 658613e5ff0Sstefano_zampini if (R_owns_row_starts) hypre_ParCSRMatrixSetRowStartsOwner(hR,1); 659613e5ff0Sstefano_zampini PetscFunctionReturn(0); 660613e5ff0Sstefano_zampini } 661613e5ff0Sstefano_zampini 662613e5ff0Sstefano_zampini #undef __FUNCT__ 6636f231fbdSstefano_zampini #define __FUNCT__ "MatPtAPNumeric_AIJ_AIJ_wHYPRE" 6646f231fbdSstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat P,Mat C) 665613e5ff0Sstefano_zampini { 6666f231fbdSstefano_zampini Mat B; 667613e5ff0Sstefano_zampini hypre_ParCSRMatrix *hA,*hP,*hPtAP; 668613e5ff0Sstefano_zampini PetscErrorCode ierr; 669613e5ff0Sstefano_zampini 670613e5ff0Sstefano_zampini PetscFunctionBegin; 671613e5ff0Sstefano_zampini ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr); 672613e5ff0Sstefano_zampini ierr = MatAIJGetParCSR_Private(P,&hP);CHKERRQ(ierr); 673613e5ff0Sstefano_zampini ierr = MatHYPRE_ParCSR_RAP(hP,hA,hP,&hPtAP);CHKERRQ(ierr); 6746f231fbdSstefano_zampini ierr = MatCreateFromParCSR(hPtAP,MATAIJ,PETSC_OWN_POINTER,&B);CHKERRQ(ierr); 6756f231fbdSstefano_zampini ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr); 676613e5ff0Sstefano_zampini ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr); 677613e5ff0Sstefano_zampini ierr = MatAIJRestoreParCSR_Private(P,&hP);CHKERRQ(ierr); 6786f231fbdSstefano_zampini PetscFunctionReturn(0); 6796f231fbdSstefano_zampini } 6806f231fbdSstefano_zampini 6816f231fbdSstefano_zampini #undef __FUNCT__ 6826f231fbdSstefano_zampini #define __FUNCT__ "MatPtAPSymbolic_AIJ_AIJ_wHYPRE" 6836f231fbdSstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat P,PetscReal fill,Mat *C) 6846f231fbdSstefano_zampini { 6856f231fbdSstefano_zampini PetscErrorCode ierr; 6866f231fbdSstefano_zampini PetscFunctionBegin; 6876f231fbdSstefano_zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr); 6886f231fbdSstefano_zampini ierr = MatSetType(*C,MATAIJ);CHKERRQ(ierr); 6896f231fbdSstefano_zampini (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_AIJ_wHYPRE; 690613e5ff0Sstefano_zampini PetscFunctionReturn(0); 691613e5ff0Sstefano_zampini } 692613e5ff0Sstefano_zampini 693613e5ff0Sstefano_zampini #undef __FUNCT__ 6944cc28894Sstefano_zampini #define __FUNCT__ "MatPtAPNumeric_AIJ_HYPRE" 6954cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_HYPRE(Mat A,Mat P,Mat C) 696613e5ff0Sstefano_zampini { 6974cc28894Sstefano_zampini Mat B; 6984cc28894Sstefano_zampini Mat_HYPRE *hP; 699613e5ff0Sstefano_zampini hypre_ParCSRMatrix *hA = NULL,*Pparcsr,*ptapparcsr; 700613e5ff0Sstefano_zampini HYPRE_Int type; 701613e5ff0Sstefano_zampini MPI_Comm comm = PetscObjectComm((PetscObject)A); 7024cc28894Sstefano_zampini PetscBool ishypre; 703613e5ff0Sstefano_zampini PetscErrorCode ierr; 704613e5ff0Sstefano_zampini 705613e5ff0Sstefano_zampini PetscFunctionBegin; 7064cc28894Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr); 7074cc28894Sstefano_zampini if (!ishypre) SETERRQ1(comm,PETSC_ERR_USER,"P should be of type %s",MATHYPRE); 7084cc28894Sstefano_zampini hP = (Mat_HYPRE*)P->data; 709613e5ff0Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type)); 710613e5ff0Sstefano_zampini if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 711613e5ff0Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr)); 712613e5ff0Sstefano_zampini 713613e5ff0Sstefano_zampini /* It looks like we don't need to have the diagonal entries 714613e5ff0Sstefano_zampini ordered first in the rows of the diagonal part 715613e5ff0Sstefano_zampini for boomerAMGBuildCoarseOperator to work */ 716613e5ff0Sstefano_zampini ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr); 717613e5ff0Sstefano_zampini ierr = MatHYPRE_ParCSR_RAP(Pparcsr,hA,Pparcsr,&ptapparcsr);CHKERRQ(ierr); 718613e5ff0Sstefano_zampini ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr); 719225daaf8SStefano Zampini 7204cc28894Sstefano_zampini /* create temporary matrix and merge to C */ 7214cc28894Sstefano_zampini ierr = MatCreateFromParCSR(ptapparcsr,((PetscObject)C)->type_name,PETSC_OWN_POINTER,&B);CHKERRQ(ierr); 7224cc28894Sstefano_zampini ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr); 7234cc28894Sstefano_zampini PetscFunctionReturn(0); 7244cc28894Sstefano_zampini } 7254cc28894Sstefano_zampini 7264cc28894Sstefano_zampini #undef __FUNCT__ 7274cc28894Sstefano_zampini #define __FUNCT__ "MatPtAP_AIJ_HYPRE" 7284cc28894Sstefano_zampini static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 7294cc28894Sstefano_zampini { 7304cc28894Sstefano_zampini PetscErrorCode ierr; 7314cc28894Sstefano_zampini 7324cc28894Sstefano_zampini PetscFunctionBegin; 7334cc28894Sstefano_zampini if (scall == MAT_INITIAL_MATRIX) { 7344cc28894Sstefano_zampini const char *deft = MATAIJ; 7354cc28894Sstefano_zampini char type[256]; 7364cc28894Sstefano_zampini PetscBool flg; 7374cc28894Sstefano_zampini 7384cc28894Sstefano_zampini ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 7394cc28894Sstefano_zampini ierr = PetscOptionsFList("-matptap_hypre_outtype","Matrix type",NULL,MatList,deft,type,256,&flg);CHKERRQ(ierr); 7404cc28894Sstefano_zampini ierr = PetscOptionsEnd();CHKERRQ(ierr); 7414cc28894Sstefano_zampini ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 7424cc28894Sstefano_zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr); 7434cc28894Sstefano_zampini if (flg) { 7444cc28894Sstefano_zampini ierr = MatSetType(*C,type);CHKERRQ(ierr); 7454cc28894Sstefano_zampini } else { 7464cc28894Sstefano_zampini ierr = MatSetType(*C,deft);CHKERRQ(ierr); 7474cc28894Sstefano_zampini } 7484cc28894Sstefano_zampini (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_HYPRE; 7494cc28894Sstefano_zampini ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 7504cc28894Sstefano_zampini } 7514cc28894Sstefano_zampini ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 7524cc28894Sstefano_zampini ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 753613e5ff0Sstefano_zampini ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 754c1a070e6SStefano Zampini PetscFunctionReturn(0); 755c1a070e6SStefano Zampini } 756c1a070e6SStefano Zampini 757c1a070e6SStefano Zampini #undef __FUNCT__ 7584cc28894Sstefano_zampini #define __FUNCT__ "MatPtAPNumeric_HYPRE_HYPRE" 7594cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_HYPRE_HYPRE(Mat A,Mat P,Mat C) 7604cc28894Sstefano_zampini { 7614cc28894Sstefano_zampini Mat B; 7624cc28894Sstefano_zampini hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr; 7634cc28894Sstefano_zampini Mat_HYPRE *hA,*hP; 7644cc28894Sstefano_zampini PetscBool ishypre; 7654cc28894Sstefano_zampini HYPRE_Int type; 7664cc28894Sstefano_zampini PetscErrorCode ierr; 7674cc28894Sstefano_zampini 7684cc28894Sstefano_zampini PetscFunctionBegin; 7694cc28894Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr); 7704cc28894Sstefano_zampini if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)P),PETSC_ERR_USER,"P should be of type %s",MATHYPRE); 7714cc28894Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr); 7724cc28894Sstefano_zampini if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE); 7734cc28894Sstefano_zampini hA = (Mat_HYPRE*)A->data; 7744cc28894Sstefano_zampini hP = (Mat_HYPRE*)P->data; 7754cc28894Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 7764cc28894Sstefano_zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 7774cc28894Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type)); 7784cc28894Sstefano_zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 7794cc28894Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr)); 7804cc28894Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr)); 7814cc28894Sstefano_zampini ierr = MatHYPRE_ParCSR_RAP(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr);CHKERRQ(ierr); 7824cc28894Sstefano_zampini ierr = MatCreateFromParCSR(ptapparcsr,MATHYPRE,PETSC_OWN_POINTER,&B);CHKERRQ(ierr); 7834cc28894Sstefano_zampini ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr); 7844cc28894Sstefano_zampini PetscFunctionReturn(0); 7854cc28894Sstefano_zampini } 7864cc28894Sstefano_zampini 7874cc28894Sstefano_zampini #undef __FUNCT__ 788cd8bc7baSStefano Zampini #define __FUNCT__ "MatPtAP_HYPRE_HYPRE" 789cd8bc7baSStefano Zampini static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 790cd8bc7baSStefano Zampini { 791cd8bc7baSStefano Zampini PetscErrorCode ierr; 792cd8bc7baSStefano Zampini 793cd8bc7baSStefano Zampini PetscFunctionBegin; 7944cc28894Sstefano_zampini if (scall == MAT_INITIAL_MATRIX) { 7954cc28894Sstefano_zampini ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 7964cc28894Sstefano_zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr); 7974cc28894Sstefano_zampini ierr = MatSetType(*C,MATHYPRE);CHKERRQ(ierr); 7984cc28894Sstefano_zampini (*C)->ops->ptapnumeric = MatPtAPNumeric_HYPRE_HYPRE; 7994cc28894Sstefano_zampini ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 8004cc28894Sstefano_zampini } 801613e5ff0Sstefano_zampini ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 8024cc28894Sstefano_zampini ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 803613e5ff0Sstefano_zampini ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 804cd8bc7baSStefano Zampini PetscFunctionReturn(0); 805cd8bc7baSStefano Zampini } 806cd8bc7baSStefano Zampini 807cd8bc7baSStefano Zampini #undef __FUNCT__ 808*5e5acdf2Sstefano_zampini #define __FUNCT__ "MatMatMultNumeric_AIJ_AIJ_wHYPRE" 809*5e5acdf2Sstefano_zampini static PetscErrorCode MatMatMultNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C) 810*5e5acdf2Sstefano_zampini { 811*5e5acdf2Sstefano_zampini Mat D; 812*5e5acdf2Sstefano_zampini hypre_ParCSRMatrix *hA,*hB,*hAB; 813*5e5acdf2Sstefano_zampini PetscErrorCode ierr; 814*5e5acdf2Sstefano_zampini 815*5e5acdf2Sstefano_zampini PetscFunctionBegin; 816*5e5acdf2Sstefano_zampini ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr); 817*5e5acdf2Sstefano_zampini ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr); 818*5e5acdf2Sstefano_zampini PetscStackPush("hypre_ParMatmul"); 819*5e5acdf2Sstefano_zampini hAB = hypre_ParMatmul(hA,hB); 820*5e5acdf2Sstefano_zampini PetscStackPop; 821*5e5acdf2Sstefano_zampini ierr = MatCreateFromParCSR(hAB,MATAIJ,PETSC_OWN_POINTER,&D);CHKERRQ(ierr); 822*5e5acdf2Sstefano_zampini ierr = MatHeaderMerge(C,&D);CHKERRQ(ierr); 823*5e5acdf2Sstefano_zampini ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr); 824*5e5acdf2Sstefano_zampini ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr); 825*5e5acdf2Sstefano_zampini PetscFunctionReturn(0); 826*5e5acdf2Sstefano_zampini } 827*5e5acdf2Sstefano_zampini 828*5e5acdf2Sstefano_zampini #undef __FUNCT__ 829*5e5acdf2Sstefano_zampini #define __FUNCT__ "MatMatMultSymbolic_AIJ_AIJ_wHYPRE" 830*5e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat B,PetscReal fill,Mat *C) 831*5e5acdf2Sstefano_zampini { 832*5e5acdf2Sstefano_zampini PetscErrorCode ierr; 833*5e5acdf2Sstefano_zampini PetscFunctionBegin; 834*5e5acdf2Sstefano_zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr); 835*5e5acdf2Sstefano_zampini ierr = MatSetType(*C,MATAIJ);CHKERRQ(ierr); 836*5e5acdf2Sstefano_zampini (*C)->ops->matmultnumeric = MatMatMultNumeric_AIJ_AIJ_wHYPRE; 837*5e5acdf2Sstefano_zampini PetscFunctionReturn(0); 838*5e5acdf2Sstefano_zampini } 839*5e5acdf2Sstefano_zampini 840*5e5acdf2Sstefano_zampini #undef __FUNCT__ 84163c07aadSStefano Zampini #define __FUNCT__ "MatMultTranspose_HYPRE" 842ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y) 84363c07aadSStefano Zampini { 84463c07aadSStefano Zampini PetscErrorCode ierr; 84563c07aadSStefano Zampini 84663c07aadSStefano Zampini PetscFunctionBegin; 84763c07aadSStefano Zampini ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr); 84863c07aadSStefano Zampini PetscFunctionReturn(0); 84963c07aadSStefano Zampini } 85063c07aadSStefano Zampini 85163c07aadSStefano Zampini #undef __FUNCT__ 85263c07aadSStefano Zampini #define __FUNCT__ "MatMult_HYPRE" 853ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y) 85463c07aadSStefano Zampini { 85563c07aadSStefano Zampini PetscErrorCode ierr; 85663c07aadSStefano Zampini 85763c07aadSStefano Zampini PetscFunctionBegin; 85863c07aadSStefano Zampini ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr); 85963c07aadSStefano Zampini PetscFunctionReturn(0); 86063c07aadSStefano Zampini } 86163c07aadSStefano Zampini 86263c07aadSStefano Zampini #undef __FUNCT__ 86363c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_MultKernel_Private" 864ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans) 86563c07aadSStefano Zampini { 86663c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 86763c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 86863c07aadSStefano Zampini hypre_ParVector *hx,*hy; 86963c07aadSStefano Zampini PetscScalar *ax,*ay,*sax,*say; 87063c07aadSStefano Zampini PetscErrorCode ierr; 87163c07aadSStefano Zampini 87263c07aadSStefano Zampini PetscFunctionBegin; 87363c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 87463c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx)); 87563c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy)); 87663c07aadSStefano Zampini ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); 87763c07aadSStefano Zampini ierr = VecGetArray(y,&ay);CHKERRQ(ierr); 87863c07aadSStefano Zampini if (trans) { 87958968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,ay,say); 88058968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax); 88163c07aadSStefano Zampini hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx); 88258968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,say,ay); 88358968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax); 88463c07aadSStefano Zampini } else { 88558968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax); 88658968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,ay,say); 88763c07aadSStefano Zampini hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy); 88858968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax); 88958968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,say,ay); 89063c07aadSStefano Zampini } 89163c07aadSStefano Zampini ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); 89263c07aadSStefano Zampini ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr); 89363c07aadSStefano Zampini PetscFunctionReturn(0); 89463c07aadSStefano Zampini } 89563c07aadSStefano Zampini 89663c07aadSStefano Zampini #undef __FUNCT__ 89763c07aadSStefano Zampini #define __FUNCT__ "MatDestroy_HYPRE" 898ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A) 89963c07aadSStefano Zampini { 90063c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 90163c07aadSStefano Zampini PetscErrorCode ierr; 90263c07aadSStefano Zampini 90363c07aadSStefano Zampini PetscFunctionBegin; 90463c07aadSStefano Zampini if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x)); 90563c07aadSStefano Zampini if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b)); 906978814f1SStefano Zampini if (hA->ij) { 907978814f1SStefano Zampini if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL; 908978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij)); 909978814f1SStefano Zampini } 91063c07aadSStefano Zampini if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); } 91163c07aadSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr); 9122df22349SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr); 913c1a070e6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr); 914c1a070e6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr); 915d975228cSstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPRESetPreallocation_C",NULL);CHKERRQ(ierr); 916dd9c0a25Sstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPREGetParCSR_C",NULL);CHKERRQ(ierr); 91763c07aadSStefano Zampini ierr = PetscFree(A->data);CHKERRQ(ierr); 91863c07aadSStefano Zampini PetscFunctionReturn(0); 91963c07aadSStefano Zampini } 92063c07aadSStefano Zampini 92163c07aadSStefano Zampini #undef __FUNCT__ 92263c07aadSStefano Zampini #define __FUNCT__ "MatSetUp_HYPRE" 923ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A) 92463c07aadSStefano Zampini { 9254ec6421dSstefano_zampini PetscErrorCode ierr; 9264ec6421dSstefano_zampini 9274ec6421dSstefano_zampini PetscFunctionBegin; 9284ec6421dSstefano_zampini ierr = MatHYPRESetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr); 9294ec6421dSstefano_zampini PetscFunctionReturn(0); 9304ec6421dSstefano_zampini } 9314ec6421dSstefano_zampini 9324ec6421dSstefano_zampini #undef __FUNCT__ 9334ec6421dSstefano_zampini #define __FUNCT__ "MatAssemblyEnd_HYPRE" 9344ec6421dSstefano_zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode) 9354ec6421dSstefano_zampini { 93663c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 93763c07aadSStefano Zampini Vec x,b; 93863c07aadSStefano Zampini PetscErrorCode ierr; 93963c07aadSStefano Zampini 94063c07aadSStefano Zampini PetscFunctionBegin; 9414ec6421dSstefano_zampini if (mode == MAT_FLUSH_ASSEMBLY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_FLUSH_ASSEMBLY currently not supported with MATHYPRE"); 9424ec6421dSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij)); 94363c07aadSStefano Zampini if (hA->x) PetscFunctionReturn(0); 94463c07aadSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 94563c07aadSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 94663c07aadSStefano Zampini ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr); 94763c07aadSStefano Zampini ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr); 94863c07aadSStefano Zampini ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr); 94963c07aadSStefano Zampini ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr); 95063c07aadSStefano Zampini ierr = VecDestroy(&x);CHKERRQ(ierr); 95163c07aadSStefano Zampini ierr = VecDestroy(&b);CHKERRQ(ierr); 95263c07aadSStefano Zampini PetscFunctionReturn(0); 95363c07aadSStefano Zampini } 95463c07aadSStefano Zampini 955d975228cSstefano_zampini #define MATHYPRE_SCRATCH 2048 956d975228cSstefano_zampini 957d975228cSstefano_zampini #undef __FUNCT__ 958d975228cSstefano_zampini #define __FUNCT__ "MatSetValues_HYPRE" 959d975228cSstefano_zampini PetscErrorCode MatSetValues_HYPRE(Mat A, PetscInt nr, const PetscInt rows[], PetscInt nc, const PetscInt cols[], const PetscScalar v[], InsertMode ins) 960d975228cSstefano_zampini { 961d975228cSstefano_zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 962d975228cSstefano_zampini PetscScalar *vals = (PetscScalar *)v; 963d975228cSstefano_zampini PetscScalar sscr[MATHYPRE_SCRATCH]; 9647d968826Sstefano_zampini HYPRE_Int cscr[2][MATHYPRE_SCRATCH]; 9657d968826Sstefano_zampini HYPRE_Int i,nzc; 966d975228cSstefano_zampini PetscErrorCode ierr; 967d975228cSstefano_zampini 968d975228cSstefano_zampini PetscFunctionBegin; 969d975228cSstefano_zampini for (i=0,nzc=0;i<nc;i++) { 970d975228cSstefano_zampini if (cols[i] >= 0) { 971d975228cSstefano_zampini cscr[0][nzc ] = cols[i]; 972d975228cSstefano_zampini cscr[1][nzc++] = i; 973d975228cSstefano_zampini } 974d975228cSstefano_zampini } 975d975228cSstefano_zampini if (!nzc) PetscFunctionReturn(0); 976d975228cSstefano_zampini 977d975228cSstefano_zampini if (ins == ADD_VALUES) { 978d975228cSstefano_zampini for (i=0;i<nr;i++) { 979d975228cSstefano_zampini if (rows[i] >= 0) { 980d975228cSstefano_zampini PetscInt j; 981d975228cSstefano_zampini for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]]; 9827d968826Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixAddToValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr)); 983d975228cSstefano_zampini } 984d975228cSstefano_zampini vals += nc; 985d975228cSstefano_zampini } 986d975228cSstefano_zampini } else { /* INSERT_VALUES */ 987d975228cSstefano_zampini #if defined(PETSC_USE_DEBUG) 988d975228cSstefano_zampini /* Insert values cannot be used to insert offproc entries */ 989d975228cSstefano_zampini PetscInt rst,ren; 990d975228cSstefano_zampini ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr); 991d975228cSstefano_zampini for (i=0;i<nr;i++) 992d975228cSstefano_zampini if (rows[i] >= 0 && (rows[i] < rst || rows[i] >= ren)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use INSERT_VALUES for off-proc entries with MatHYPRE. Use ADD_VALUES instead"); 993d975228cSstefano_zampini #endif 994d975228cSstefano_zampini for (i=0;i<nr;i++) { 995d975228cSstefano_zampini if (rows[i] >= 0) { 996d975228cSstefano_zampini PetscInt j; 997d975228cSstefano_zampini for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]]; 9987d968826Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixSetValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr)); 999d975228cSstefano_zampini } 1000d975228cSstefano_zampini vals += nc; 1001d975228cSstefano_zampini } 1002d975228cSstefano_zampini } 1003d975228cSstefano_zampini PetscFunctionReturn(0); 1004d975228cSstefano_zampini } 1005d975228cSstefano_zampini 1006d975228cSstefano_zampini #undef __FUNCT__ 1007d975228cSstefano_zampini #define __FUNCT__ "MatHYPRESetPreallocation_HYPRE" 1008d975228cSstefano_zampini static PetscErrorCode MatHYPRESetPreallocation_HYPRE(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[]) 1009d975228cSstefano_zampini { 1010d975228cSstefano_zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 10117d968826Sstefano_zampini HYPRE_Int *hdnnz,*honnz; 101206a29025Sstefano_zampini PetscInt i,rs,re,cs,ce,bs; 1013d975228cSstefano_zampini PetscMPIInt size; 1014d975228cSstefano_zampini PetscErrorCode ierr; 1015d975228cSstefano_zampini 1016d975228cSstefano_zampini PetscFunctionBegin; 101706a29025Sstefano_zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1018d975228cSstefano_zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1019d975228cSstefano_zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1020d975228cSstefano_zampini rs = A->rmap->rstart; 1021d975228cSstefano_zampini re = A->rmap->rend; 1022d975228cSstefano_zampini cs = A->cmap->rstart; 1023d975228cSstefano_zampini ce = A->cmap->rend; 1024d975228cSstefano_zampini if (!hA->ij) { 1025d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rs,re-1,cs,ce-1,&hA->ij)); 1026d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); 1027d975228cSstefano_zampini } else { 10287d968826Sstefano_zampini HYPRE_Int hrs,hre,hcs,hce; 1029d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetLocalRange,(hA->ij,&hrs,&hre,&hcs,&hce)); 1030d975228cSstefano_zampini if (hre-hrs+1 != re -rs) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent local rows: IJMatrix [%D,%D), PETSc [%D,%d)",hrs,hre+1,rs,re); 1031d975228cSstefano_zampini if (hce-hcs+1 != ce -cs) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent local cols: IJMatrix [%D,%D), PETSc [%D,%d)",hcs,hce+1,cs,ce); 1032d975228cSstefano_zampini } 1033d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij)); 1034d975228cSstefano_zampini 103506a29025Sstefano_zampini if (dnz == PETSC_DEFAULT || dnz == PETSC_DECIDE) dnz = 10*bs; 103606a29025Sstefano_zampini if (onz == PETSC_DEFAULT || onz == PETSC_DECIDE) onz = 10*bs; 103706a29025Sstefano_zampini 1038d975228cSstefano_zampini if (!dnnz) { 1039d975228cSstefano_zampini ierr = PetscMalloc1(A->rmap->n,&hdnnz);CHKERRQ(ierr); 1040d975228cSstefano_zampini for (i=0;i<A->rmap->n;i++) hdnnz[i] = dnz; 1041d975228cSstefano_zampini } else { 10427d968826Sstefano_zampini hdnnz = (HYPRE_Int*)dnnz; 1043d975228cSstefano_zampini } 1044d975228cSstefano_zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 1045d975228cSstefano_zampini if (size > 1) { 1046d975228cSstefano_zampini if (!onnz) { 1047d975228cSstefano_zampini ierr = PetscMalloc1(A->rmap->n,&honnz);CHKERRQ(ierr); 1048d975228cSstefano_zampini for (i=0;i<A->rmap->n;i++) honnz[i] = onz; 1049d975228cSstefano_zampini } else { 10507d968826Sstefano_zampini honnz = (HYPRE_Int*)onnz; 1051d975228cSstefano_zampini } 1052d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(hA->ij,hdnnz,honnz)); 1053d975228cSstefano_zampini } else { 1054d975228cSstefano_zampini honnz = NULL; 1055d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixSetRowSizes,(hA->ij,hdnnz)); 1056d975228cSstefano_zampini } 1057d975228cSstefano_zampini if (!dnnz) { 1058d975228cSstefano_zampini ierr = PetscFree(hdnnz);CHKERRQ(ierr); 1059d975228cSstefano_zampini } 1060d975228cSstefano_zampini if (!onnz && honnz) { 1061d975228cSstefano_zampini ierr = PetscFree(honnz);CHKERRQ(ierr); 1062d975228cSstefano_zampini } 106306a29025Sstefano_zampini A->preallocated = PETSC_TRUE; 1064d975228cSstefano_zampini 1065d975228cSstefano_zampini /* SetDiagOffdSizes sets hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0 */ 1066d975228cSstefano_zampini { 1067d975228cSstefano_zampini hypre_AuxParCSRMatrix *aux_matrix; 1068d975228cSstefano_zampini aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij); 1069d975228cSstefano_zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1; 1070d975228cSstefano_zampini } 1071d975228cSstefano_zampini PetscFunctionReturn(0); 1072d975228cSstefano_zampini } 1073d975228cSstefano_zampini 1074d975228cSstefano_zampini /*@C 1075d975228cSstefano_zampini MatHYPRESetPreallocation - Preallocates memory for a sparse parallel matrix in HYPRE IJ format 1076d975228cSstefano_zampini 1077d975228cSstefano_zampini Collective on Mat 1078d975228cSstefano_zampini 1079d975228cSstefano_zampini Input Parameters: 1080d975228cSstefano_zampini + A - the matrix 1081d975228cSstefano_zampini . dnz - number of nonzeros per row in DIAGONAL portion of local submatrix 1082d975228cSstefano_zampini (same value is used for all local rows) 1083d975228cSstefano_zampini . dnnz - array containing the number of nonzeros in the various rows of the 1084d975228cSstefano_zampini DIAGONAL portion of the local submatrix (possibly different for each row) 1085d975228cSstefano_zampini or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure. 1086d975228cSstefano_zampini The size of this array is equal to the number of local rows, i.e 'm'. 1087d975228cSstefano_zampini For matrices that will be factored, you must leave room for (and set) 1088d975228cSstefano_zampini the diagonal entry even if it is zero. 1089d975228cSstefano_zampini . onz - number of nonzeros per row in the OFF-DIAGONAL portion of local 1090d975228cSstefano_zampini submatrix (same value is used for all local rows). 1091d975228cSstefano_zampini - onnz - array containing the number of nonzeros in the various rows of the 1092d975228cSstefano_zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 1093d975228cSstefano_zampini each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero 1094d975228cSstefano_zampini structure. The size of this array is equal to the number 1095d975228cSstefano_zampini of local rows, i.e 'm'. 1096d975228cSstefano_zampini 1097d975228cSstefano_zampini Notes: If the *nnz parameter is given then the *nz parameter is ignored; for sequential matrices, onz and onnz are ignored. 1098d975228cSstefano_zampini 1099d975228cSstefano_zampini Level: intermediate 1100d975228cSstefano_zampini 1101d975228cSstefano_zampini .keywords: matrix, aij, compressed row, sparse, parallel 1102d975228cSstefano_zampini 1103d975228cSstefano_zampini .seealso: MatCreate(), MatMPIAIJSetPreallocation, MATHYPRE 1104d975228cSstefano_zampini @*/ 1105d975228cSstefano_zampini #undef __FUNCT__ 1106d975228cSstefano_zampini #define __FUNCT__ "MatHYPRESetPreallocation" 1107d975228cSstefano_zampini PetscErrorCode MatHYPRESetPreallocation(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[]) 1108d975228cSstefano_zampini { 1109d975228cSstefano_zampini PetscErrorCode ierr; 1110d975228cSstefano_zampini 1111d975228cSstefano_zampini PetscFunctionBegin; 1112d975228cSstefano_zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1113d975228cSstefano_zampini PetscValidType(A,1); 1114d975228cSstefano_zampini ierr = PetscTryMethod(A,"MatHYPRESetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(A,dnz,dnnz,onz,onnz));CHKERRQ(ierr); 1115d975228cSstefano_zampini PetscFunctionReturn(0); 1116d975228cSstefano_zampini } 1117d975228cSstefano_zampini 1118225daaf8SStefano Zampini /* 1119225daaf8SStefano Zampini MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix 1120225daaf8SStefano Zampini 1121225daaf8SStefano Zampini Collective 1122225daaf8SStefano Zampini 1123225daaf8SStefano Zampini Input Parameters: 1124225daaf8SStefano Zampini + vparcsr - the pointer to the hypre_ParCSRMatrix 1125bb4689ddSStefano Zampini . mtype - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported. 1126225daaf8SStefano Zampini - copymode - PETSc copying options 1127225daaf8SStefano Zampini 1128225daaf8SStefano Zampini Output Parameter: 1129225daaf8SStefano Zampini . A - the matrix 1130225daaf8SStefano Zampini 1131225daaf8SStefano Zampini Level: intermediate 1132225daaf8SStefano Zampini 1133225daaf8SStefano Zampini .seealso: MatHYPRE, PetscCopyMode 1134225daaf8SStefano Zampini */ 1135978814f1SStefano Zampini #undef __FUNCT__ 1136225daaf8SStefano Zampini #define __FUNCT__ "MatCreateFromParCSR" 1137dd9c0a25Sstefano_zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(hypre_ParCSRMatrix *vparcsr, MatType mtype, PetscCopyMode copymode, Mat* A) 1138978814f1SStefano Zampini { 1139225daaf8SStefano Zampini Mat T; 1140978814f1SStefano Zampini Mat_HYPRE *hA; 114163c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 1142978814f1SStefano Zampini MPI_Comm comm; 1143978814f1SStefano Zampini PetscInt rstart,rend,cstart,cend,M,N; 1144bb4689ddSStefano Zampini PetscBool isseqaij,ismpiaij,isaij,ishyp,isis; 1145978814f1SStefano Zampini PetscErrorCode ierr; 1146978814f1SStefano Zampini 1147978814f1SStefano Zampini PetscFunctionBegin; 1148978814f1SStefano Zampini parcsr = (hypre_ParCSRMatrix *)vparcsr; 1149978814f1SStefano Zampini comm = hypre_ParCSRMatrixComm(parcsr); 1150225daaf8SStefano Zampini ierr = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 1151225daaf8SStefano Zampini ierr = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 1152225daaf8SStefano Zampini ierr = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr); 1153225daaf8SStefano Zampini ierr = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr); 11548cfe8d00SStefano Zampini ierr = PetscStrcmp(mtype,MATIS,&isis);CHKERRQ(ierr); 1155225daaf8SStefano Zampini isaij = (PetscBool)(isseqaij || ismpiaij || isaij); 1156bb4689ddSStefano Zampini if (!isaij && !ishyp && !isis) SETERRQ6(comm,PETSC_ERR_SUP,"Unsupported MatType %s! Supported types are %s, %s, %s, %s, and %s",mtype,MATAIJ,MATSEQAIJ,MATMPIAIJ,MATIS,MATHYPRE); 1157225daaf8SStefano Zampini if (ishyp && copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES"); 1158978814f1SStefano Zampini 1159978814f1SStefano Zampini /* access ParCSRMatrix */ 1160978814f1SStefano Zampini rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr); 1161978814f1SStefano Zampini rend = hypre_ParCSRMatrixLastRowIndex(parcsr); 1162978814f1SStefano Zampini cstart = hypre_ParCSRMatrixFirstColDiag(parcsr); 1163978814f1SStefano Zampini cend = hypre_ParCSRMatrixLastColDiag(parcsr); 1164978814f1SStefano Zampini M = hypre_ParCSRMatrixGlobalNumRows(parcsr); 1165978814f1SStefano Zampini N = hypre_ParCSRMatrixGlobalNumCols(parcsr); 1166978814f1SStefano Zampini 1167fa92c42cSstefano_zampini /* fix for empty local rows/columns */ 1168fa92c42cSstefano_zampini if (rend < rstart) rend = rstart; 1169fa92c42cSstefano_zampini if (cend < cstart) cend = cstart; 1170fa92c42cSstefano_zampini 1171978814f1SStefano Zampini /* create PETSc matrix with MatHYPRE */ 1172225daaf8SStefano Zampini ierr = MatCreate(comm,&T);CHKERRQ(ierr); 1173225daaf8SStefano Zampini ierr = MatSetSizes(T,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr); 1174225daaf8SStefano Zampini ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr); 1175225daaf8SStefano Zampini hA = (Mat_HYPRE*)(T->data); 1176978814f1SStefano Zampini 1177978814f1SStefano Zampini /* create HYPRE_IJMatrix */ 1178978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij)); 1179978814f1SStefano Zampini 1180978814f1SStefano Zampini /* set ParCSR object */ 1181978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); 1182978814f1SStefano Zampini hypre_IJMatrixObject(hA->ij) = parcsr; 11834ec6421dSstefano_zampini T->preallocated = PETSC_TRUE; 1184978814f1SStefano Zampini 1185978814f1SStefano Zampini /* set assembled flag */ 1186978814f1SStefano Zampini hypre_IJMatrixAssembleFlag(hA->ij) = 1; 1187978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij)); 1188225daaf8SStefano Zampini if (ishyp) { 11896d2a658fSstefano_zampini PetscMPIInt myid = 0; 11906d2a658fSstefano_zampini 11916d2a658fSstefano_zampini /* make sure we always have row_starts and col_starts available */ 11926d2a658fSstefano_zampini if (HYPRE_AssumedPartitionCheck()) { 11936d2a658fSstefano_zampini ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr); 11946d2a658fSstefano_zampini } 11956d2a658fSstefano_zampini if (!hypre_ParCSRMatrixOwnsColStarts(parcsr)) { 11966d2a658fSstefano_zampini PetscLayout map; 11976d2a658fSstefano_zampini 11986d2a658fSstefano_zampini ierr = MatGetLayouts(T,NULL,&map);CHKERRQ(ierr); 11996d2a658fSstefano_zampini ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 12006d2a658fSstefano_zampini hypre_ParCSRMatrixColStarts(parcsr) = map->range + myid; 12016d2a658fSstefano_zampini } 12026d2a658fSstefano_zampini if (!hypre_ParCSRMatrixOwnsRowStarts(parcsr)) { 12036d2a658fSstefano_zampini PetscLayout map; 12046d2a658fSstefano_zampini 12056d2a658fSstefano_zampini ierr = MatGetLayouts(T,&map,NULL);CHKERRQ(ierr); 12066d2a658fSstefano_zampini ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 12076d2a658fSstefano_zampini hypre_ParCSRMatrixRowStarts(parcsr) = map->range + myid; 12086d2a658fSstefano_zampini } 1209978814f1SStefano Zampini /* prevent from freeing the pointer */ 1210978814f1SStefano Zampini if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE; 1211225daaf8SStefano Zampini *A = T; 1212978814f1SStefano Zampini ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1213978814f1SStefano Zampini ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1214bb4689ddSStefano Zampini } else if (isaij) { 1215bb4689ddSStefano Zampini if (copymode != PETSC_OWN_POINTER) { 1216225daaf8SStefano Zampini /* prevent from freeing the pointer */ 1217225daaf8SStefano Zampini hA->inner_free = PETSC_FALSE; 1218225daaf8SStefano Zampini ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr); 1219225daaf8SStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 1220225daaf8SStefano Zampini } else { /* AIJ return type with PETSC_OWN_POINTER */ 1221225daaf8SStefano Zampini ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr); 1222225daaf8SStefano Zampini *A = T; 1223225daaf8SStefano Zampini } 1224bb4689ddSStefano Zampini } else if (isis) { 12258cfe8d00SStefano Zampini ierr = MatConvert_HYPRE_IS(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr); 12268cfe8d00SStefano Zampini if (copymode != PETSC_OWN_POINTER) hA->inner_free = PETSC_FALSE; 1227bb4689ddSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 1228bb4689ddSStefano Zampini } 1229978814f1SStefano Zampini PetscFunctionReturn(0); 1230978814f1SStefano Zampini } 1231978814f1SStefano Zampini 123263c07aadSStefano Zampini #undef __FUNCT__ 1233dd9c0a25Sstefano_zampini #define __FUNCT__ "MatHYPREGetParCSR_HYPRE" 1234dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR_HYPRE(Mat A, hypre_ParCSRMatrix **parcsr) 1235dd9c0a25Sstefano_zampini { 1236dd9c0a25Sstefano_zampini Mat_HYPRE* hA = (Mat_HYPRE*)A->data; 1237dd9c0a25Sstefano_zampini HYPRE_Int type; 1238dd9c0a25Sstefano_zampini PetscErrorCode ierr; 1239dd9c0a25Sstefano_zampini 1240dd9c0a25Sstefano_zampini PetscFunctionBegin; 1241dd9c0a25Sstefano_zampini if (!hA->ij) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"HYPRE_IJMatrix not present"); 1242dd9c0a25Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 1243dd9c0a25Sstefano_zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"HYPRE_IJMatrix is not of type HYPRE_PARCSR"); 1244dd9c0a25Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)parcsr)); 1245dd9c0a25Sstefano_zampini PetscFunctionReturn(0); 1246dd9c0a25Sstefano_zampini } 1247dd9c0a25Sstefano_zampini 1248dd9c0a25Sstefano_zampini /* 1249dd9c0a25Sstefano_zampini MatHYPREGetParCSR - Gets the pointer to the ParCSR matrix 1250dd9c0a25Sstefano_zampini 1251dd9c0a25Sstefano_zampini Not collective 1252dd9c0a25Sstefano_zampini 1253dd9c0a25Sstefano_zampini Input Parameters: 1254dd9c0a25Sstefano_zampini + A - the MATHYPRE object 1255dd9c0a25Sstefano_zampini 1256dd9c0a25Sstefano_zampini Output Parameter: 1257dd9c0a25Sstefano_zampini . parcsr - the pointer to the hypre_ParCSRMatrix 1258dd9c0a25Sstefano_zampini 1259dd9c0a25Sstefano_zampini Level: intermediate 1260dd9c0a25Sstefano_zampini 1261dd9c0a25Sstefano_zampini .seealso: MatHYPRE, PetscCopyMode 1262dd9c0a25Sstefano_zampini */ 1263dd9c0a25Sstefano_zampini #undef __FUNCT__ 1264dd9c0a25Sstefano_zampini #define __FUNCT__ "MatHYPREGetParCSR" 1265dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR(Mat A, hypre_ParCSRMatrix **parcsr) 1266dd9c0a25Sstefano_zampini { 1267dd9c0a25Sstefano_zampini PetscErrorCode ierr; 1268dd9c0a25Sstefano_zampini 1269dd9c0a25Sstefano_zampini PetscFunctionBegin; 1270dd9c0a25Sstefano_zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1271dd9c0a25Sstefano_zampini PetscValidType(A,1); 1272dd9c0a25Sstefano_zampini ierr = PetscUseMethod(A,"MatHYPREGetParCSR_C",(Mat,hypre_ParCSRMatrix**),(A,parcsr));CHKERRQ(ierr); 1273dd9c0a25Sstefano_zampini PetscFunctionReturn(0); 1274dd9c0a25Sstefano_zampini } 1275dd9c0a25Sstefano_zampini 1276dd9c0a25Sstefano_zampini #undef __FUNCT__ 127763c07aadSStefano Zampini #define __FUNCT__ "MatCreate_HYPRE" 127863c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B) 127963c07aadSStefano Zampini { 128063c07aadSStefano Zampini Mat_HYPRE *hB; 128163c07aadSStefano Zampini PetscErrorCode ierr; 128263c07aadSStefano Zampini 128363c07aadSStefano Zampini PetscFunctionBegin; 128463c07aadSStefano Zampini ierr = PetscNewLog(B,&hB);CHKERRQ(ierr); 1285978814f1SStefano Zampini hB->inner_free = PETSC_TRUE; 1286978814f1SStefano Zampini 128763c07aadSStefano Zampini B->data = (void*)hB; 128863c07aadSStefano Zampini B->rmap->bs = 1; 128963c07aadSStefano Zampini B->assembled = PETSC_FALSE; 129063c07aadSStefano Zampini 129163c07aadSStefano Zampini B->ops->mult = MatMult_HYPRE; 129263c07aadSStefano Zampini B->ops->multtranspose = MatMultTranspose_HYPRE; 129363c07aadSStefano Zampini B->ops->setup = MatSetUp_HYPRE; 129463c07aadSStefano Zampini B->ops->destroy = MatDestroy_HYPRE; 129563c07aadSStefano Zampini B->ops->assemblyend = MatAssemblyEnd_HYPRE; 1296cd8bc7baSStefano Zampini B->ops->ptap = MatPtAP_HYPRE_HYPRE; 1297d975228cSstefano_zampini B->ops->setvalues = MatSetValues_HYPRE; 129863c07aadSStefano Zampini 129963c07aadSStefano Zampini ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr); 130063c07aadSStefano Zampini ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr); 130163c07aadSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr); 13022df22349SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr); 1303c1a070e6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr); 1304c1a070e6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr); 1305d975228cSstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPRESetPreallocation_C",MatHYPRESetPreallocation_HYPRE);CHKERRQ(ierr); 1306dd9c0a25Sstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPREGetParCSR_C",MatHYPREGetParCSR_HYPRE);CHKERRQ(ierr); 130763c07aadSStefano Zampini PetscFunctionReturn(0); 130863c07aadSStefano Zampini } 130963c07aadSStefano Zampini 1310225daaf8SStefano Zampini #undef __FUNCT__ 1311225daaf8SStefano Zampini #define __FUNCT__ "hypre_array_destroy" 1312225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void *ptr) 1313225daaf8SStefano Zampini { 1314225daaf8SStefano Zampini PetscFunctionBegin; 1315225daaf8SStefano Zampini hypre_TFree(ptr); 1316225daaf8SStefano Zampini PetscFunctionReturn(0); 1317225daaf8SStefano Zampini } 1318225daaf8SStefano Zampini 131963c07aadSStefano Zampini #if 0 132063c07aadSStefano Zampini /* 132163c07aadSStefano Zampini Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format 132263c07aadSStefano Zampini 132363c07aadSStefano Zampini This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first 132463c07aadSStefano Zampini which will corrupt the PETSc data structure if we did this. Need a work around to this problem. 132563c07aadSStefano Zampini */ 132663c07aadSStefano Zampini #include <_hypre_IJ_mv.h> 132763c07aadSStefano Zampini #include <HYPRE_IJ_mv.h> 132863c07aadSStefano Zampini 132963c07aadSStefano Zampini #undef __FUNCT__ 133063c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixLink" 133163c07aadSStefano Zampini PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij) 133263c07aadSStefano Zampini { 133363c07aadSStefano Zampini PetscErrorCode ierr; 133463c07aadSStefano Zampini PetscInt rstart,rend,cstart,cend; 133563c07aadSStefano Zampini PetscBool flg; 133663c07aadSStefano Zampini hypre_AuxParCSRMatrix *aux_matrix; 133763c07aadSStefano Zampini 133863c07aadSStefano Zampini PetscFunctionBegin; 133963c07aadSStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 134063c07aadSStefano Zampini PetscValidType(A,1); 134163c07aadSStefano Zampini PetscValidPointer(ij,2); 134263c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 134363c07aadSStefano Zampini if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices"); 134463c07aadSStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 134563c07aadSStefano Zampini 134663c07aadSStefano Zampini ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 134763c07aadSStefano Zampini rstart = A->rmap->rstart; 134863c07aadSStefano Zampini rend = A->rmap->rend; 134963c07aadSStefano Zampini cstart = A->cmap->rstart; 135063c07aadSStefano Zampini cend = A->cmap->rend; 135163c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij)); 135263c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR)); 135363c07aadSStefano Zampini 135463c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij)); 135563c07aadSStefano Zampini PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij)); 135663c07aadSStefano Zampini 135763c07aadSStefano Zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 135863c07aadSStefano Zampini 135963c07aadSStefano Zampini /* this is the Hack part where we monkey directly with the hypre datastructures */ 136063c07aadSStefano Zampini 136163c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij)); 136263c07aadSStefano Zampini ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 136363c07aadSStefano Zampini PetscFunctionReturn(0); 136463c07aadSStefano Zampini } 136563c07aadSStefano Zampini #endif 1366