163c07aadSStefano Zampini 263c07aadSStefano Zampini /* 363c07aadSStefano Zampini Creates hypre ijmatrix from PETSc matrix 463c07aadSStefano Zampini */ 563c07aadSStefano Zampini #include <petscsys.h> 663c07aadSStefano Zampini #include <petsc/private/matimpl.h> 763c07aadSStefano Zampini #include <../src/mat/impls/hypre/mhypre.h> 863c07aadSStefano Zampini #include <../src/mat/impls/aij/mpi/mpiaij.h> 958968eb6SStefano Zampini #include <../src/vec/vec/impls/hypre/vhyp.h> 1058968eb6SStefano Zampini #include <HYPRE.h> 11*cd8bc7baSStefano Zampini #include <_hypre_parcsr_ls.h> 1263c07aadSStefano Zampini 1363c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*); 1463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix); 1563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix); 1663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix); 1763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool); 1863c07aadSStefano Zampini 1963c07aadSStefano Zampini #undef __FUNCT__ 2063c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate" 2163c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij) 2263c07aadSStefano Zampini { 2363c07aadSStefano Zampini PetscErrorCode ierr; 2463c07aadSStefano Zampini PetscInt i,n_d,n_o; 2563c07aadSStefano Zampini const PetscInt *ia_d,*ia_o; 2663c07aadSStefano Zampini PetscBool done_d=PETSC_FALSE,done_o=PETSC_FALSE; 2763c07aadSStefano Zampini PetscInt *nnz_d=NULL,*nnz_o=NULL; 2863c07aadSStefano Zampini 2963c07aadSStefano Zampini PetscFunctionBegin; 3063c07aadSStefano Zampini if (A_d) { /* determine number of nonzero entries in local diagonal part */ 3163c07aadSStefano Zampini ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr); 3263c07aadSStefano Zampini if (done_d) { 3363c07aadSStefano Zampini ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr); 3463c07aadSStefano Zampini for (i=0; i<n_d; i++) { 3563c07aadSStefano Zampini nnz_d[i] = ia_d[i+1] - ia_d[i]; 3663c07aadSStefano Zampini } 3763c07aadSStefano Zampini } 3863c07aadSStefano Zampini ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr); 3963c07aadSStefano Zampini } 4063c07aadSStefano Zampini if (A_o) { /* determine number of nonzero entries in local off-diagonal part */ 4163c07aadSStefano Zampini ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr); 4263c07aadSStefano Zampini if (done_o) { 4363c07aadSStefano Zampini ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr); 4463c07aadSStefano Zampini for (i=0; i<n_o; i++) { 4563c07aadSStefano Zampini nnz_o[i] = ia_o[i+1] - ia_o[i]; 4663c07aadSStefano Zampini } 4763c07aadSStefano Zampini } 4863c07aadSStefano Zampini ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr); 4963c07aadSStefano Zampini } 5063c07aadSStefano Zampini if (done_d) { /* set number of nonzeros in HYPRE IJ matrix */ 5163c07aadSStefano Zampini if (!done_o) { /* only diagonal part */ 5263c07aadSStefano Zampini ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr); 5363c07aadSStefano Zampini for (i=0; i<n_d; i++) { 5463c07aadSStefano Zampini nnz_o[i] = 0; 5563c07aadSStefano Zampini } 5663c07aadSStefano Zampini } 5763c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o)); 5863c07aadSStefano Zampini ierr = PetscFree(nnz_d);CHKERRQ(ierr); 5963c07aadSStefano Zampini ierr = PetscFree(nnz_o);CHKERRQ(ierr); 6063c07aadSStefano Zampini } 6163c07aadSStefano Zampini PetscFunctionReturn(0); 6263c07aadSStefano Zampini } 6363c07aadSStefano Zampini 6463c07aadSStefano Zampini #undef __FUNCT__ 6563c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_CreateFromMat" 6663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA) 6763c07aadSStefano Zampini { 6863c07aadSStefano Zampini PetscErrorCode ierr; 6963c07aadSStefano Zampini PetscInt rstart,rend,cstart,cend; 7063c07aadSStefano Zampini 7163c07aadSStefano Zampini PetscFunctionBegin; 7263c07aadSStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 7363c07aadSStefano Zampini rstart = A->rmap->rstart; 7463c07aadSStefano Zampini rend = A->rmap->rend; 7563c07aadSStefano Zampini cstart = A->cmap->rstart; 7663c07aadSStefano Zampini cend = A->cmap->rend; 7763c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij)); 7863c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); 7963c07aadSStefano Zampini { 8063c07aadSStefano Zampini PetscBool same; 8163c07aadSStefano Zampini Mat A_d,A_o; 8263c07aadSStefano Zampini const PetscInt *colmap; 8363c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr); 8463c07aadSStefano Zampini if (same) { 8563c07aadSStefano Zampini ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); 8663c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr); 8763c07aadSStefano Zampini PetscFunctionReturn(0); 8863c07aadSStefano Zampini } 8963c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr); 9063c07aadSStefano Zampini if (same) { 9163c07aadSStefano Zampini ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); 9263c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr); 9363c07aadSStefano Zampini PetscFunctionReturn(0); 9463c07aadSStefano Zampini } 9563c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr); 9663c07aadSStefano Zampini if (same) { 9763c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr); 9863c07aadSStefano Zampini PetscFunctionReturn(0); 9963c07aadSStefano Zampini } 10063c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr); 10163c07aadSStefano Zampini if (same) { 10263c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr); 10363c07aadSStefano Zampini PetscFunctionReturn(0); 10463c07aadSStefano Zampini } 10563c07aadSStefano Zampini } 10663c07aadSStefano Zampini PetscFunctionReturn(0); 10763c07aadSStefano Zampini } 10863c07aadSStefano Zampini 10963c07aadSStefano Zampini #undef __FUNCT__ 11063c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixCopy" 11163c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij) 11263c07aadSStefano Zampini { 11363c07aadSStefano Zampini PetscErrorCode ierr; 11463c07aadSStefano Zampini PetscInt i,rstart,rend,ncols,nr,nc; 11563c07aadSStefano Zampini const PetscScalar *values; 11663c07aadSStefano Zampini const PetscInt *cols; 11763c07aadSStefano Zampini PetscBool flg; 11863c07aadSStefano Zampini 11963c07aadSStefano Zampini PetscFunctionBegin; 12063c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 12163c07aadSStefano Zampini ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr); 12263c07aadSStefano Zampini if (flg && nr == nc) { 12363c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr); 12463c07aadSStefano Zampini PetscFunctionReturn(0); 12563c07aadSStefano Zampini } 12663c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); 12763c07aadSStefano Zampini if (flg) { 12863c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr); 12963c07aadSStefano Zampini PetscFunctionReturn(0); 13063c07aadSStefano Zampini } 13163c07aadSStefano Zampini 13263c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 13363c07aadSStefano Zampini ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 13463c07aadSStefano Zampini for (i=rstart; i<rend; i++) { 13563c07aadSStefano Zampini ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 13663c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values)); 13763c07aadSStefano Zampini ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 13863c07aadSStefano Zampini } 13963c07aadSStefano Zampini PetscFunctionReturn(0); 14063c07aadSStefano Zampini } 14163c07aadSStefano Zampini 14263c07aadSStefano Zampini #undef __FUNCT__ 14363c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ" 14463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij) 14563c07aadSStefano Zampini { 14663c07aadSStefano Zampini PetscErrorCode ierr; 14763c07aadSStefano Zampini Mat_SeqAIJ *pdiag = (Mat_SeqAIJ*)A->data; 14858968eb6SStefano Zampini HYPRE_Int type; 14963c07aadSStefano Zampini hypre_ParCSRMatrix *par_matrix; 15063c07aadSStefano Zampini hypre_AuxParCSRMatrix *aux_matrix; 15163c07aadSStefano Zampini hypre_CSRMatrix *hdiag; 15263c07aadSStefano Zampini 15363c07aadSStefano Zampini PetscFunctionBegin; 15463c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 155ea9daf28SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type)); 156ea9daf28SStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 157ea9daf28SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix)); 15863c07aadSStefano Zampini hdiag = hypre_ParCSRMatrixDiag(par_matrix); 15963c07aadSStefano Zampini /* 16063c07aadSStefano Zampini this is the Hack part where we monkey directly with the hypre datastructures 16163c07aadSStefano Zampini */ 16263c07aadSStefano Zampini ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt)); 16363c07aadSStefano Zampini ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt)); 16463c07aadSStefano Zampini ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); 165ea9daf28SStefano Zampini 166ea9daf28SStefano Zampini aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); 16763c07aadSStefano Zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 16863c07aadSStefano Zampini PetscFunctionReturn(0); 16963c07aadSStefano Zampini } 17063c07aadSStefano Zampini 17163c07aadSStefano Zampini #undef __FUNCT__ 17263c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ" 17363c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij) 17463c07aadSStefano Zampini { 17563c07aadSStefano Zampini PetscErrorCode ierr; 17663c07aadSStefano Zampini Mat_MPIAIJ *pA = (Mat_MPIAIJ*)A->data; 17763c07aadSStefano Zampini Mat_SeqAIJ *pdiag,*poffd; 17863c07aadSStefano Zampini PetscInt i,*garray = pA->garray,*jj,cstart,*pjj; 17958968eb6SStefano Zampini HYPRE_Int type; 18063c07aadSStefano Zampini hypre_ParCSRMatrix *par_matrix; 18163c07aadSStefano Zampini hypre_AuxParCSRMatrix *aux_matrix; 18263c07aadSStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 18363c07aadSStefano Zampini 18463c07aadSStefano Zampini PetscFunctionBegin; 18563c07aadSStefano Zampini pdiag = (Mat_SeqAIJ*) pA->A->data; 18663c07aadSStefano Zampini poffd = (Mat_SeqAIJ*) pA->B->data; 18763c07aadSStefano Zampini /* cstart is only valid for square MPIAIJ layed out in the usual way */ 18863c07aadSStefano Zampini ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr); 18963c07aadSStefano Zampini 19063c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 191ea9daf28SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type)); 192ea9daf28SStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 193ea9daf28SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix)); 19463c07aadSStefano Zampini hdiag = hypre_ParCSRMatrixDiag(par_matrix); 19563c07aadSStefano Zampini hoffd = hypre_ParCSRMatrixOffd(par_matrix); 19663c07aadSStefano Zampini 19763c07aadSStefano Zampini /* 19863c07aadSStefano Zampini this is the Hack part where we monkey directly with the hypre datastructures 19963c07aadSStefano Zampini */ 20063c07aadSStefano Zampini ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); 20163c07aadSStefano Zampini /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */ 20263c07aadSStefano Zampini jj = (PetscInt*)hdiag->j; 20363c07aadSStefano Zampini pjj = (PetscInt*)pdiag->j; 20463c07aadSStefano Zampini for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i]; 20563c07aadSStefano Zampini ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar)); 20663c07aadSStefano Zampini ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt)); 20763c07aadSStefano Zampini /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this 20863c07aadSStefano Zampini If we hacked a hypre a bit more we might be able to avoid this step */ 20963c07aadSStefano Zampini jj = (PetscInt*) hoffd->j; 21063c07aadSStefano Zampini pjj = (PetscInt*) poffd->j; 21163c07aadSStefano Zampini for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]]; 21263c07aadSStefano Zampini ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar)); 21363c07aadSStefano Zampini 214ea9daf28SStefano Zampini aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); 21563c07aadSStefano Zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 21663c07aadSStefano Zampini PetscFunctionReturn(0); 21763c07aadSStefano Zampini } 21863c07aadSStefano Zampini 21963c07aadSStefano Zampini #undef __FUNCT__ 2202df22349SStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_IS" 2212df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B) 2222df22349SStefano Zampini { 2232df22349SStefano Zampini Mat_HYPRE* mhA = (Mat_HYPRE*)(A->data); 2242df22349SStefano Zampini Mat lA; 2252df22349SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 2262df22349SStefano Zampini IS is; 2272df22349SStefano Zampini hypre_ParCSRMatrix *hA; 2282df22349SStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 2292df22349SStefano Zampini MPI_Comm comm; 2302df22349SStefano Zampini PetscScalar *hdd,*hod,*aa,*data; 2312df22349SStefano Zampini PetscInt *col_map_offd,*hdi,*hdj,*hoi,*hoj; 2322df22349SStefano Zampini PetscInt *ii,*jj,*iptr,*jptr; 2332df22349SStefano Zampini PetscInt cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N; 23458968eb6SStefano Zampini HYPRE_Int type; 2352df22349SStefano Zampini PetscErrorCode ierr; 2362df22349SStefano Zampini 2372df22349SStefano Zampini PetscFunctionBegin; 238a1787963SStefano Zampini comm = PetscObjectComm((PetscObject)A); 2392df22349SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type)); 2402df22349SStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 2412df22349SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA)); 2422df22349SStefano Zampini M = hypre_ParCSRMatrixGlobalNumRows(hA); 2432df22349SStefano Zampini N = hypre_ParCSRMatrixGlobalNumCols(hA); 2442df22349SStefano Zampini str = hypre_ParCSRMatrixFirstRowIndex(hA); 2452df22349SStefano Zampini stc = hypre_ParCSRMatrixFirstColDiag(hA); 2462df22349SStefano Zampini hdiag = hypre_ParCSRMatrixDiag(hA); 2472df22349SStefano Zampini hoffd = hypre_ParCSRMatrixOffd(hA); 2482df22349SStefano Zampini dr = hypre_CSRMatrixNumRows(hdiag); 2492df22349SStefano Zampini dc = hypre_CSRMatrixNumCols(hdiag); 2502df22349SStefano Zampini nnz = hypre_CSRMatrixNumNonzeros(hdiag); 2512df22349SStefano Zampini hdi = hypre_CSRMatrixI(hdiag); 2522df22349SStefano Zampini hdj = hypre_CSRMatrixJ(hdiag); 2532df22349SStefano Zampini hdd = hypre_CSRMatrixData(hdiag); 2542df22349SStefano Zampini oc = hypre_CSRMatrixNumCols(hoffd); 2552df22349SStefano Zampini nnz += hypre_CSRMatrixNumNonzeros(hoffd); 2562df22349SStefano Zampini hoi = hypre_CSRMatrixI(hoffd); 2572df22349SStefano Zampini hoj = hypre_CSRMatrixJ(hoffd); 2582df22349SStefano Zampini hod = hypre_CSRMatrixData(hoffd); 2592df22349SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 2602df22349SStefano Zampini PetscInt *aux; 2612df22349SStefano Zampini 2622df22349SStefano Zampini /* generate l2g maps for rows and cols */ 2632df22349SStefano Zampini ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr); 2642df22349SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 2652df22349SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 2662df22349SStefano Zampini col_map_offd = hypre_ParCSRMatrixColMapOffd(hA); 2672df22349SStefano Zampini ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr); 2682df22349SStefano Zampini for (i=0; i<dc; i++) aux[i] = i+stc; 2692df22349SStefano Zampini for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i]; 2702df22349SStefano Zampini ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 2712df22349SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 2722df22349SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 2732df22349SStefano Zampini /* create MATIS object */ 2742df22349SStefano Zampini ierr = MatCreate(comm,B);CHKERRQ(ierr); 2752df22349SStefano Zampini ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr); 2762df22349SStefano Zampini ierr = MatSetType(*B,MATIS);CHKERRQ(ierr); 2772df22349SStefano Zampini ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr); 2782df22349SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 2792df22349SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 2802df22349SStefano Zampini 2812df22349SStefano Zampini /* allocate CSR for local matrix */ 2822df22349SStefano Zampini ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr); 2832df22349SStefano Zampini ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr); 2842df22349SStefano Zampini ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr); 2852df22349SStefano Zampini } else { 2862df22349SStefano Zampini PetscInt nr; 2872df22349SStefano Zampini PetscBool done; 2882df22349SStefano Zampini ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr); 2892df22349SStefano Zampini ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr); 2902df22349SStefano 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); 2912df22349SStefano 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); 2922df22349SStefano Zampini ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr); 2932df22349SStefano Zampini } 2942df22349SStefano Zampini /* merge local matrices */ 2952df22349SStefano Zampini ii = iptr; 2962df22349SStefano Zampini jj = jptr; 2972df22349SStefano Zampini aa = data; 2982df22349SStefano Zampini *ii = *(hdi++) + *(hoi++); 2992df22349SStefano Zampini for (jd=0,jo=0,cum=0; *ii<nnz; cum++) { 3002df22349SStefano Zampini PetscScalar *aold = aa; 3012df22349SStefano Zampini PetscInt *jold = jj,nc = jd+jo; 3022df22349SStefano Zampini for (; jd<*hdi; jd++) { *jj++ = *hdj++; *aa++ = *hdd++; } 3032df22349SStefano Zampini for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; } 3042df22349SStefano Zampini *(++ii) = *(hdi++) + *(hoi++); 3052df22349SStefano Zampini ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr); 3062df22349SStefano Zampini } 3072df22349SStefano Zampini for (; cum<dr; cum++) *(++ii) = nnz; 3082df22349SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 309a033916dSStefano Zampini Mat_SeqAIJ* a; 310a033916dSStefano Zampini 3112df22349SStefano Zampini ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr); 3122df22349SStefano Zampini ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr); 313a033916dSStefano Zampini /* hack SeqAIJ */ 314a033916dSStefano Zampini a = (Mat_SeqAIJ*)(lA->data); 315a033916dSStefano Zampini a->free_a = PETSC_TRUE; 316a033916dSStefano Zampini a->free_ij = PETSC_TRUE; 3172df22349SStefano Zampini ierr = MatDestroy(&lA);CHKERRQ(ierr); 3182df22349SStefano Zampini } 3192df22349SStefano Zampini ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3202df22349SStefano Zampini ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3212df22349SStefano Zampini if (reuse == MAT_INPLACE_MATRIX) { 3222df22349SStefano Zampini ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 3232df22349SStefano Zampini } 3242df22349SStefano Zampini PetscFunctionReturn(0); 3252df22349SStefano Zampini } 3262df22349SStefano Zampini 3272df22349SStefano Zampini #undef __FUNCT__ 32863c07aadSStefano Zampini #define __FUNCT__ "MatConvert_AIJ_HYPRE" 32963c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B) 33063c07aadSStefano Zampini { 33163c07aadSStefano Zampini Mat_HYPRE *hB; 33263c07aadSStefano Zampini MPI_Comm comm = PetscObjectComm((PetscObject)A); 33363c07aadSStefano Zampini PetscErrorCode ierr; 33463c07aadSStefano Zampini 33563c07aadSStefano Zampini PetscFunctionBegin; 33663c07aadSStefano Zampini if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX"); 33763c07aadSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 33863c07aadSStefano Zampini /* always destroy the old matrix and create a new memory; 33963c07aadSStefano Zampini hope this does not churn the memory too much. The problem 34063c07aadSStefano Zampini is I do not know if it is possible to put the matrix back to 34163c07aadSStefano Zampini its initial state so that we can directly copy the values 34263c07aadSStefano Zampini the second time through. */ 34363c07aadSStefano Zampini hB = (Mat_HYPRE*)((*B)->data); 34463c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij)); 34563c07aadSStefano Zampini } else { 34663c07aadSStefano Zampini ierr = MatCreate(comm,B);CHKERRQ(ierr); 34763c07aadSStefano Zampini ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr); 34863c07aadSStefano Zampini ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 34963c07aadSStefano Zampini ierr = MatSetUp(*B);CHKERRQ(ierr); 35063c07aadSStefano Zampini hB = (Mat_HYPRE*)((*B)->data); 35163c07aadSStefano Zampini } 35263c07aadSStefano Zampini ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr); 35363c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr); 35463c07aadSStefano Zampini ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 35563c07aadSStefano Zampini ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 35663c07aadSStefano Zampini PetscFunctionReturn(0); 35763c07aadSStefano Zampini } 35863c07aadSStefano Zampini 35963c07aadSStefano Zampini #undef __FUNCT__ 36063c07aadSStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_AIJ" 361ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B) 36263c07aadSStefano Zampini { 36363c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 36463c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 36563c07aadSStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 36663c07aadSStefano Zampini MPI_Comm comm; 36763c07aadSStefano Zampini PetscScalar *da,*oa,*aptr; 36863c07aadSStefano Zampini PetscInt *dii,*djj,*oii,*ojj,*iptr; 36963c07aadSStefano Zampini PetscInt i,dnnz,onnz,m,n; 37058968eb6SStefano Zampini HYPRE_Int type; 37163c07aadSStefano Zampini PetscMPIInt size; 37263c07aadSStefano Zampini PetscErrorCode ierr; 37363c07aadSStefano Zampini 37463c07aadSStefano Zampini PetscFunctionBegin; 37563c07aadSStefano Zampini comm = PetscObjectComm((PetscObject)A); 37663c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 37763c07aadSStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 37863c07aadSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 37963c07aadSStefano Zampini PetscBool ismpiaij,isseqaij; 38063c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 38163c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 38263c07aadSStefano Zampini if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported"); 38363c07aadSStefano Zampini } 38463c07aadSStefano Zampini ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 38563c07aadSStefano Zampini 38663c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 38763c07aadSStefano Zampini hdiag = hypre_ParCSRMatrixDiag(parcsr); 38863c07aadSStefano Zampini hoffd = hypre_ParCSRMatrixOffd(parcsr); 38963c07aadSStefano Zampini m = hypre_CSRMatrixNumRows(hdiag); 39063c07aadSStefano Zampini n = hypre_CSRMatrixNumCols(hdiag); 39163c07aadSStefano Zampini dnnz = hypre_CSRMatrixNumNonzeros(hdiag); 39263c07aadSStefano Zampini onnz = hypre_CSRMatrixNumNonzeros(hoffd); 39363c07aadSStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 39463c07aadSStefano Zampini ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr); 39563c07aadSStefano Zampini ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr); 39663c07aadSStefano Zampini ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr); 39763c07aadSStefano Zampini } else { 39863c07aadSStefano Zampini PetscInt nr; 39963c07aadSStefano Zampini PetscBool done; 40063c07aadSStefano Zampini if (size > 1) { 40163c07aadSStefano Zampini Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); 40263c07aadSStefano Zampini 40363c07aadSStefano Zampini ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); 40463c07aadSStefano 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); 40563c07aadSStefano 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); 40663c07aadSStefano Zampini ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr); 40763c07aadSStefano Zampini } else { 40863c07aadSStefano Zampini ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); 40963c07aadSStefano Zampini if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m); 41063c07aadSStefano 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); 41163c07aadSStefano Zampini ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr); 41263c07aadSStefano Zampini } 41363c07aadSStefano Zampini } 41463c07aadSStefano Zampini ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 41563c07aadSStefano Zampini ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr); 41663c07aadSStefano Zampini ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr); 41763c07aadSStefano Zampini iptr = djj; 41863c07aadSStefano Zampini aptr = da; 41963c07aadSStefano Zampini for (i=0; i<m; i++) { 42063c07aadSStefano Zampini PetscInt nc = dii[i+1]-dii[i]; 42163c07aadSStefano Zampini ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr); 42263c07aadSStefano Zampini iptr += nc; 42363c07aadSStefano Zampini aptr += nc; 42463c07aadSStefano Zampini } 42563c07aadSStefano Zampini if (size > 1) { 42663c07aadSStefano Zampini PetscInt *offdj,*coffd; 42763c07aadSStefano Zampini 42863c07aadSStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 42963c07aadSStefano Zampini ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr); 43063c07aadSStefano Zampini ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr); 43163c07aadSStefano Zampini ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr); 43263c07aadSStefano Zampini } else { 43363c07aadSStefano Zampini Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); 43463c07aadSStefano Zampini PetscInt nr,hr = hypre_CSRMatrixNumRows(hoffd); 43563c07aadSStefano Zampini PetscBool done; 43663c07aadSStefano Zampini 43763c07aadSStefano Zampini ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr); 43863c07aadSStefano 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); 43963c07aadSStefano 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); 44063c07aadSStefano Zampini ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr); 44163c07aadSStefano Zampini } 44263c07aadSStefano Zampini ierr = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 44363c07aadSStefano Zampini offdj = hypre_CSRMatrixJ(hoffd); 44463c07aadSStefano Zampini coffd = hypre_ParCSRMatrixColMapOffd(parcsr); 44563c07aadSStefano Zampini for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]]; 44663c07aadSStefano Zampini ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr); 44763c07aadSStefano Zampini iptr = ojj; 44863c07aadSStefano Zampini aptr = oa; 44963c07aadSStefano Zampini for (i=0; i<m; i++) { 45063c07aadSStefano Zampini PetscInt nc = oii[i+1]-oii[i]; 45163c07aadSStefano Zampini ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr); 45263c07aadSStefano Zampini iptr += nc; 45363c07aadSStefano Zampini aptr += nc; 45463c07aadSStefano Zampini } 45563c07aadSStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 45663c07aadSStefano Zampini Mat_MPIAIJ *b; 45763c07aadSStefano Zampini Mat_SeqAIJ *d,*o; 45863c07aadSStefano Zampini ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii, 45963c07aadSStefano Zampini djj,da,oii,ojj,oa,B);CHKERRQ(ierr); 46063c07aadSStefano Zampini /* hack MPIAIJ */ 46163c07aadSStefano Zampini b = (Mat_MPIAIJ*)((*B)->data); 46263c07aadSStefano Zampini d = (Mat_SeqAIJ*)b->A->data; 46363c07aadSStefano Zampini o = (Mat_SeqAIJ*)b->B->data; 46463c07aadSStefano Zampini d->free_a = PETSC_TRUE; 46563c07aadSStefano Zampini d->free_ij = PETSC_TRUE; 46663c07aadSStefano Zampini o->free_a = PETSC_TRUE; 46763c07aadSStefano Zampini o->free_ij = PETSC_TRUE; 46863c07aadSStefano Zampini } 46963c07aadSStefano Zampini } else if (reuse != MAT_REUSE_MATRIX) { 47063c07aadSStefano Zampini Mat_SeqAIJ* b; 47163c07aadSStefano Zampini ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr); 47263c07aadSStefano Zampini /* hack SeqAIJ */ 47363c07aadSStefano Zampini b = (Mat_SeqAIJ*)((*B)->data); 47463c07aadSStefano Zampini b->free_a = PETSC_TRUE; 47563c07aadSStefano Zampini b->free_ij = PETSC_TRUE; 47663c07aadSStefano Zampini } 47763c07aadSStefano Zampini if (reuse == MAT_INPLACE_MATRIX) { 47863c07aadSStefano Zampini ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 47963c07aadSStefano Zampini } 48063c07aadSStefano Zampini PetscFunctionReturn(0); 48163c07aadSStefano Zampini } 48263c07aadSStefano Zampini 48363c07aadSStefano Zampini #undef __FUNCT__ 484*cd8bc7baSStefano Zampini #define __FUNCT__ "MatPtAP_HYPRE_HYPRE" 485*cd8bc7baSStefano Zampini static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 486*cd8bc7baSStefano Zampini { 487*cd8bc7baSStefano Zampini hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr; 488*cd8bc7baSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data, *hP = (Mat_HYPRE*)P->data; 489*cd8bc7baSStefano Zampini HYPRE_Int type; 490*cd8bc7baSStefano Zampini PetscErrorCode ierr; 491*cd8bc7baSStefano Zampini 492*cd8bc7baSStefano Zampini PetscFunctionBegin; 493*cd8bc7baSStefano Zampini if (scall == MAT_REUSE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Unsupported MAT_REUSE_MATRIX"); 494*cd8bc7baSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 495*cd8bc7baSStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 496*cd8bc7baSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type)); 497*cd8bc7baSStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 498*cd8bc7baSStefano Zampini ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 499*cd8bc7baSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr)); 500*cd8bc7baSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr)); 501*cd8bc7baSStefano Zampini PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr)); 502*cd8bc7baSStefano Zampini ierr = MatHYPRECreateFromParCSR(ptapparcsr,PETSC_USE_POINTER,C);CHKERRQ(ierr); 503*cd8bc7baSStefano Zampini ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 504*cd8bc7baSStefano Zampini PetscFunctionReturn(0); 505*cd8bc7baSStefano Zampini } 506*cd8bc7baSStefano Zampini 507*cd8bc7baSStefano Zampini #undef __FUNCT__ 50863c07aadSStefano Zampini #define __FUNCT__ "MatMultTranspose_HYPRE" 509ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y) 51063c07aadSStefano Zampini { 51163c07aadSStefano Zampini PetscErrorCode ierr; 51263c07aadSStefano Zampini 51363c07aadSStefano Zampini PetscFunctionBegin; 51463c07aadSStefano Zampini ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr); 51563c07aadSStefano Zampini PetscFunctionReturn(0); 51663c07aadSStefano Zampini } 51763c07aadSStefano Zampini 51863c07aadSStefano Zampini #undef __FUNCT__ 51963c07aadSStefano Zampini #define __FUNCT__ "MatMult_HYPRE" 520ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y) 52163c07aadSStefano Zampini { 52263c07aadSStefano Zampini PetscErrorCode ierr; 52363c07aadSStefano Zampini 52463c07aadSStefano Zampini PetscFunctionBegin; 52563c07aadSStefano Zampini ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr); 52663c07aadSStefano Zampini PetscFunctionReturn(0); 52763c07aadSStefano Zampini } 52863c07aadSStefano Zampini 52963c07aadSStefano Zampini #undef __FUNCT__ 53063c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_MultKernel_Private" 531ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans) 53263c07aadSStefano Zampini { 53363c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 53463c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 53563c07aadSStefano Zampini hypre_ParVector *hx,*hy; 53663c07aadSStefano Zampini PetscScalar *ax,*ay,*sax,*say; 53763c07aadSStefano Zampini PetscErrorCode ierr; 53863c07aadSStefano Zampini 53963c07aadSStefano Zampini PetscFunctionBegin; 54063c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 54163c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx)); 54263c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy)); 54363c07aadSStefano Zampini ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); 54463c07aadSStefano Zampini ierr = VecGetArray(y,&ay);CHKERRQ(ierr); 54563c07aadSStefano Zampini if (trans) { 54658968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,ay,say); 54758968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax); 54863c07aadSStefano Zampini hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx); 54958968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,say,ay); 55058968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax); 55163c07aadSStefano Zampini } else { 55258968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax); 55358968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,ay,say); 55463c07aadSStefano Zampini hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy); 55558968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax); 55658968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,say,ay); 55763c07aadSStefano Zampini } 55863c07aadSStefano Zampini ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); 55963c07aadSStefano Zampini ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr); 56063c07aadSStefano Zampini PetscFunctionReturn(0); 56163c07aadSStefano Zampini } 56263c07aadSStefano Zampini 56363c07aadSStefano Zampini #undef __FUNCT__ 56463c07aadSStefano Zampini #define __FUNCT__ "MatDestroy_HYPRE" 565ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A) 56663c07aadSStefano Zampini { 56763c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 56863c07aadSStefano Zampini PetscErrorCode ierr; 56963c07aadSStefano Zampini 57063c07aadSStefano Zampini PetscFunctionBegin; 57163c07aadSStefano Zampini if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x)); 57263c07aadSStefano Zampini if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b)); 573978814f1SStefano Zampini if (hA->ij) { 574978814f1SStefano Zampini if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL; 575978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij)); 576978814f1SStefano Zampini } 57763c07aadSStefano Zampini if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); } 57863c07aadSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr); 5792df22349SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr); 58063c07aadSStefano Zampini ierr = PetscFree(A->data);CHKERRQ(ierr); 58163c07aadSStefano Zampini PetscFunctionReturn(0); 58263c07aadSStefano Zampini } 58363c07aadSStefano Zampini 58463c07aadSStefano Zampini #undef __FUNCT__ 58563c07aadSStefano Zampini #define __FUNCT__ "MatSetUp_HYPRE" 586ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A) 58763c07aadSStefano Zampini { 58863c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 58963c07aadSStefano Zampini Vec x,b; 59063c07aadSStefano Zampini PetscErrorCode ierr; 59163c07aadSStefano Zampini 59263c07aadSStefano Zampini PetscFunctionBegin; 59363c07aadSStefano Zampini if (hA->x) PetscFunctionReturn(0); 59463c07aadSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 59563c07aadSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 59663c07aadSStefano Zampini ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr); 59763c07aadSStefano Zampini ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr); 59863c07aadSStefano Zampini ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr); 59963c07aadSStefano Zampini ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr); 60063c07aadSStefano Zampini ierr = VecDestroy(&x);CHKERRQ(ierr); 60163c07aadSStefano Zampini ierr = VecDestroy(&b);CHKERRQ(ierr); 60263c07aadSStefano Zampini PetscFunctionReturn(0); 60363c07aadSStefano Zampini } 60463c07aadSStefano Zampini 60563c07aadSStefano Zampini #undef __FUNCT__ 60663c07aadSStefano Zampini #define __FUNCT__ "MatAssemblyEnd_HYPRE" 607ea9daf28SStefano Zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode) 60863c07aadSStefano Zampini { 60963c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 61063c07aadSStefano Zampini PetscErrorCode ierr; 61163c07aadSStefano Zampini 61263c07aadSStefano Zampini PetscFunctionBegin; 61363c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij)); 61463c07aadSStefano Zampini PetscFunctionReturn(0); 61563c07aadSStefano Zampini } 61663c07aadSStefano Zampini 617978814f1SStefano Zampini #undef __FUNCT__ 618978814f1SStefano Zampini #define __FUNCT__ "MatHYPRECreateFromParCSR" 619978814f1SStefano Zampini PETSC_EXTERN PetscErrorCode MatHYPRECreateFromParCSR(void *vparcsr, PetscCopyMode copymode, Mat* A) 620978814f1SStefano Zampini { 621978814f1SStefano Zampini Mat_HYPRE *hA; 62263c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 623978814f1SStefano Zampini MPI_Comm comm; 624978814f1SStefano Zampini PetscInt rstart,rend,cstart,cend,M,N; 625978814f1SStefano Zampini PetscErrorCode ierr; 626978814f1SStefano Zampini 627978814f1SStefano Zampini PetscFunctionBegin; 628978814f1SStefano Zampini parcsr = (hypre_ParCSRMatrix *)vparcsr; 629978814f1SStefano Zampini comm = hypre_ParCSRMatrixComm(parcsr); 630978814f1SStefano Zampini if (copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES"); 631978814f1SStefano Zampini 632978814f1SStefano Zampini /* access ParCSRMatrix */ 633978814f1SStefano Zampini rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr); 634978814f1SStefano Zampini rend = hypre_ParCSRMatrixLastRowIndex(parcsr); 635978814f1SStefano Zampini cstart = hypre_ParCSRMatrixFirstColDiag(parcsr); 636978814f1SStefano Zampini cend = hypre_ParCSRMatrixLastColDiag(parcsr); 637978814f1SStefano Zampini M = hypre_ParCSRMatrixGlobalNumRows(parcsr); 638978814f1SStefano Zampini N = hypre_ParCSRMatrixGlobalNumCols(parcsr); 639978814f1SStefano Zampini 640978814f1SStefano Zampini /* create PETSc matrix with MatHYPRE */ 641978814f1SStefano Zampini ierr = MatCreate(comm,A);CHKERRQ(ierr); 642978814f1SStefano Zampini ierr = MatSetSizes(*A,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr); 643978814f1SStefano Zampini ierr = MatSetType(*A,MATHYPRE);CHKERRQ(ierr); 644978814f1SStefano Zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 645978814f1SStefano Zampini hA = (Mat_HYPRE*)((*A)->data); 646978814f1SStefano Zampini 647978814f1SStefano Zampini /* create HYPRE_IJMatrix */ 648978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij)); 649978814f1SStefano Zampini 650978814f1SStefano Zampini /* set ParCSR object */ 651978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); 652978814f1SStefano Zampini hypre_IJMatrixObject(hA->ij) = parcsr; 653978814f1SStefano Zampini 654978814f1SStefano Zampini /* set assembled flag */ 655978814f1SStefano Zampini hypre_IJMatrixAssembleFlag(hA->ij) = 1; 656978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij)); 657978814f1SStefano Zampini 658978814f1SStefano Zampini /* prevent from freeing the pointer */ 659978814f1SStefano Zampini if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE; 660978814f1SStefano Zampini 661978814f1SStefano Zampini ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 662978814f1SStefano Zampini ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 663978814f1SStefano Zampini PetscFunctionReturn(0); 664978814f1SStefano Zampini } 665978814f1SStefano Zampini 66663c07aadSStefano Zampini #undef __FUNCT__ 66763c07aadSStefano Zampini #define __FUNCT__ "MatCreate_HYPRE" 66863c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B) 66963c07aadSStefano Zampini { 67063c07aadSStefano Zampini Mat_HYPRE *hB; 67163c07aadSStefano Zampini PetscErrorCode ierr; 67263c07aadSStefano Zampini 67363c07aadSStefano Zampini PetscFunctionBegin; 67463c07aadSStefano Zampini ierr = PetscNewLog(B,&hB);CHKERRQ(ierr); 675978814f1SStefano Zampini hB->inner_free = PETSC_TRUE; 676978814f1SStefano Zampini 67763c07aadSStefano Zampini B->data = (void*)hB; 67863c07aadSStefano Zampini B->rmap->bs = 1; 67963c07aadSStefano Zampini B->assembled = PETSC_FALSE; 68063c07aadSStefano Zampini 68163c07aadSStefano Zampini B->ops->mult = MatMult_HYPRE; 68263c07aadSStefano Zampini B->ops->multtranspose = MatMultTranspose_HYPRE; 68363c07aadSStefano Zampini B->ops->setup = MatSetUp_HYPRE; 68463c07aadSStefano Zampini B->ops->destroy = MatDestroy_HYPRE; 68563c07aadSStefano Zampini B->ops->assemblyend = MatAssemblyEnd_HYPRE; 686*cd8bc7baSStefano Zampini B->ops->ptap = MatPtAP_HYPRE_HYPRE; 68763c07aadSStefano Zampini 68863c07aadSStefano Zampini ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr); 68963c07aadSStefano Zampini ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr); 69063c07aadSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr); 6912df22349SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr); 69263c07aadSStefano Zampini PetscFunctionReturn(0); 69363c07aadSStefano Zampini } 69463c07aadSStefano Zampini 69563c07aadSStefano Zampini #if 0 69663c07aadSStefano Zampini /* 69763c07aadSStefano Zampini Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format 69863c07aadSStefano Zampini 69963c07aadSStefano Zampini This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first 70063c07aadSStefano Zampini which will corrupt the PETSc data structure if we did this. Need a work around to this problem. 70163c07aadSStefano Zampini */ 70263c07aadSStefano Zampini #include <_hypre_IJ_mv.h> 70363c07aadSStefano Zampini #include <HYPRE_IJ_mv.h> 70463c07aadSStefano Zampini 70563c07aadSStefano Zampini #undef __FUNCT__ 70663c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixLink" 70763c07aadSStefano Zampini PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij) 70863c07aadSStefano Zampini { 70963c07aadSStefano Zampini PetscErrorCode ierr; 71063c07aadSStefano Zampini PetscInt rstart,rend,cstart,cend; 71163c07aadSStefano Zampini PetscBool flg; 71263c07aadSStefano Zampini hypre_AuxParCSRMatrix *aux_matrix; 71363c07aadSStefano Zampini 71463c07aadSStefano Zampini PetscFunctionBegin; 71563c07aadSStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 71663c07aadSStefano Zampini PetscValidType(A,1); 71763c07aadSStefano Zampini PetscValidPointer(ij,2); 71863c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 71963c07aadSStefano Zampini if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices"); 72063c07aadSStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 72163c07aadSStefano Zampini 72263c07aadSStefano Zampini ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 72363c07aadSStefano Zampini rstart = A->rmap->rstart; 72463c07aadSStefano Zampini rend = A->rmap->rend; 72563c07aadSStefano Zampini cstart = A->cmap->rstart; 72663c07aadSStefano Zampini cend = A->cmap->rend; 72763c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij)); 72863c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR)); 72963c07aadSStefano Zampini 73063c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij)); 73163c07aadSStefano Zampini PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij)); 73263c07aadSStefano Zampini 73363c07aadSStefano Zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 73463c07aadSStefano Zampini 73563c07aadSStefano Zampini /* this is the Hack part where we monkey directly with the hypre datastructures */ 73663c07aadSStefano Zampini 73763c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij)); 73863c07aadSStefano Zampini ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 73963c07aadSStefano Zampini PetscFunctionReturn(0); 74063c07aadSStefano Zampini } 74163c07aadSStefano Zampini #endif 742