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> 963c07aadSStefano Zampini #include <HYPRE_struct_mv.h> 1063c07aadSStefano Zampini #include <HYPRE_struct_ls.h> 1163c07aadSStefano Zampini #include <_hypre_struct_mv.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; 148ea9daf28SStefano Zampini 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; 179ea9daf28SStefano Zampini 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__ 220*2df22349SStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_IS" 221*2df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B) 222*2df22349SStefano Zampini { 223*2df22349SStefano Zampini Mat_HYPRE* mhA = (Mat_HYPRE*)(A->data); 224*2df22349SStefano Zampini Mat lA; 225*2df22349SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 226*2df22349SStefano Zampini IS is; 227*2df22349SStefano Zampini hypre_ParCSRMatrix *hA; 228*2df22349SStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 229*2df22349SStefano Zampini MPI_Comm comm; 230*2df22349SStefano Zampini PetscScalar *hdd,*hod,*aa,*data; 231*2df22349SStefano Zampini PetscInt *col_map_offd,*hdi,*hdj,*hoi,*hoj; 232*2df22349SStefano Zampini PetscInt *ii,*jj,*iptr,*jptr; 233*2df22349SStefano Zampini PetscInt cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N; 234*2df22349SStefano Zampini int type; 235*2df22349SStefano Zampini PetscErrorCode ierr; 236*2df22349SStefano Zampini 237*2df22349SStefano Zampini PetscFunctionBegin; 238*2df22349SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type)); 239*2df22349SStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 240*2df22349SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA)); 241*2df22349SStefano Zampini M = hypre_ParCSRMatrixGlobalNumRows(hA); 242*2df22349SStefano Zampini N = hypre_ParCSRMatrixGlobalNumCols(hA); 243*2df22349SStefano Zampini str = hypre_ParCSRMatrixFirstRowIndex(hA); 244*2df22349SStefano Zampini stc = hypre_ParCSRMatrixFirstColDiag(hA); 245*2df22349SStefano Zampini hdiag = hypre_ParCSRMatrixDiag(hA); 246*2df22349SStefano Zampini hoffd = hypre_ParCSRMatrixOffd(hA); 247*2df22349SStefano Zampini dr = hypre_CSRMatrixNumRows(hdiag); 248*2df22349SStefano Zampini dc = hypre_CSRMatrixNumCols(hdiag); 249*2df22349SStefano Zampini nnz = hypre_CSRMatrixNumNonzeros(hdiag); 250*2df22349SStefano Zampini hdi = hypre_CSRMatrixI(hdiag); 251*2df22349SStefano Zampini hdj = hypre_CSRMatrixJ(hdiag); 252*2df22349SStefano Zampini hdd = hypre_CSRMatrixData(hdiag); 253*2df22349SStefano Zampini oc = hypre_CSRMatrixNumCols(hoffd); 254*2df22349SStefano Zampini nnz += hypre_CSRMatrixNumNonzeros(hoffd); 255*2df22349SStefano Zampini hoi = hypre_CSRMatrixI(hoffd); 256*2df22349SStefano Zampini hoj = hypre_CSRMatrixJ(hoffd); 257*2df22349SStefano Zampini hod = hypre_CSRMatrixData(hoffd); 258*2df22349SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 259*2df22349SStefano Zampini PetscInt *aux; 260*2df22349SStefano Zampini 261*2df22349SStefano Zampini /* generate l2g maps for rows and cols */ 262*2df22349SStefano Zampini comm = PetscObjectComm((PetscObject)A); 263*2df22349SStefano Zampini ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr); 264*2df22349SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 265*2df22349SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 266*2df22349SStefano Zampini col_map_offd = hypre_ParCSRMatrixColMapOffd(hA); 267*2df22349SStefano Zampini ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr); 268*2df22349SStefano Zampini for (i=0; i<dc; i++) aux[i] = i+stc; 269*2df22349SStefano Zampini for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i]; 270*2df22349SStefano Zampini ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 271*2df22349SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 272*2df22349SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 273*2df22349SStefano Zampini /* create MATIS object */ 274*2df22349SStefano Zampini ierr = MatCreate(comm,B);CHKERRQ(ierr); 275*2df22349SStefano Zampini ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr); 276*2df22349SStefano Zampini ierr = MatSetType(*B,MATIS);CHKERRQ(ierr); 277*2df22349SStefano Zampini ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr); 278*2df22349SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 279*2df22349SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 280*2df22349SStefano Zampini 281*2df22349SStefano Zampini /* allocate CSR for local matrix */ 282*2df22349SStefano Zampini ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr); 283*2df22349SStefano Zampini ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr); 284*2df22349SStefano Zampini ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr); 285*2df22349SStefano Zampini } else { 286*2df22349SStefano Zampini PetscInt nr; 287*2df22349SStefano Zampini PetscBool done; 288*2df22349SStefano Zampini ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr); 289*2df22349SStefano Zampini ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr); 290*2df22349SStefano 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); 291*2df22349SStefano 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); 292*2df22349SStefano Zampini ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr); 293*2df22349SStefano Zampini } 294*2df22349SStefano Zampini /* merge local matrices */ 295*2df22349SStefano Zampini ii = iptr; 296*2df22349SStefano Zampini jj = jptr; 297*2df22349SStefano Zampini aa = data; 298*2df22349SStefano Zampini *ii = *(hdi++) + *(hoi++); 299*2df22349SStefano Zampini for (jd=0,jo=0,cum=0; *ii<nnz; cum++) { 300*2df22349SStefano Zampini PetscScalar *aold = aa; 301*2df22349SStefano Zampini PetscInt *jold = jj,nc = jd+jo; 302*2df22349SStefano Zampini for (; jd<*hdi; jd++) { *jj++ = *hdj++; *aa++ = *hdd++; } 303*2df22349SStefano Zampini for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; } 304*2df22349SStefano Zampini *(++ii) = *(hdi++) + *(hoi++); 305*2df22349SStefano Zampini ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr); 306*2df22349SStefano Zampini } 307*2df22349SStefano Zampini for (; cum<dr; cum++) *(++ii) = nnz; 308*2df22349SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 309*2df22349SStefano Zampini ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr); 310*2df22349SStefano Zampini ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr); 311*2df22349SStefano Zampini ierr = MatDestroy(&lA);CHKERRQ(ierr); 312*2df22349SStefano Zampini } 313*2df22349SStefano Zampini ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 314*2df22349SStefano Zampini ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 315*2df22349SStefano Zampini if (reuse == MAT_INPLACE_MATRIX) { 316*2df22349SStefano Zampini ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 317*2df22349SStefano Zampini } 318*2df22349SStefano Zampini PetscFunctionReturn(0); 319*2df22349SStefano Zampini } 320*2df22349SStefano Zampini 321*2df22349SStefano Zampini #undef __FUNCT__ 32263c07aadSStefano Zampini #define __FUNCT__ "MatConvert_AIJ_HYPRE" 32363c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B) 32463c07aadSStefano Zampini { 32563c07aadSStefano Zampini Mat_HYPRE *hB; 32663c07aadSStefano Zampini MPI_Comm comm = PetscObjectComm((PetscObject)A); 32763c07aadSStefano Zampini PetscErrorCode ierr; 32863c07aadSStefano Zampini 32963c07aadSStefano Zampini PetscFunctionBegin; 33063c07aadSStefano Zampini if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX"); 33163c07aadSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 33263c07aadSStefano Zampini /* always destroy the old matrix and create a new memory; 33363c07aadSStefano Zampini hope this does not churn the memory too much. The problem 33463c07aadSStefano Zampini is I do not know if it is possible to put the matrix back to 33563c07aadSStefano Zampini its initial state so that we can directly copy the values 33663c07aadSStefano Zampini the second time through. */ 33763c07aadSStefano Zampini hB = (Mat_HYPRE*)((*B)->data); 33863c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij)); 33963c07aadSStefano Zampini } else { 34063c07aadSStefano Zampini ierr = MatCreate(comm,B);CHKERRQ(ierr); 34163c07aadSStefano Zampini ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr); 34263c07aadSStefano Zampini ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 34363c07aadSStefano Zampini ierr = MatSetUp(*B);CHKERRQ(ierr); 34463c07aadSStefano Zampini hB = (Mat_HYPRE*)((*B)->data); 34563c07aadSStefano Zampini } 34663c07aadSStefano Zampini ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr); 34763c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr); 34863c07aadSStefano Zampini ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 34963c07aadSStefano Zampini ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 35063c07aadSStefano Zampini PetscFunctionReturn(0); 35163c07aadSStefano Zampini } 35263c07aadSStefano Zampini 35363c07aadSStefano Zampini #undef __FUNCT__ 35463c07aadSStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_AIJ" 355ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B) 35663c07aadSStefano Zampini { 35763c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 35863c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 35963c07aadSStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 36063c07aadSStefano Zampini MPI_Comm comm; 36163c07aadSStefano Zampini PetscScalar *da,*oa,*aptr; 36263c07aadSStefano Zampini PetscInt *dii,*djj,*oii,*ojj,*iptr; 36363c07aadSStefano Zampini PetscInt i,dnnz,onnz,m,n; 36463c07aadSStefano Zampini int type; 36563c07aadSStefano Zampini PetscMPIInt size; 36663c07aadSStefano Zampini PetscErrorCode ierr; 36763c07aadSStefano Zampini 36863c07aadSStefano Zampini PetscFunctionBegin; 36963c07aadSStefano Zampini comm = PetscObjectComm((PetscObject)A); 37063c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 37163c07aadSStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 37263c07aadSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 37363c07aadSStefano Zampini PetscBool ismpiaij,isseqaij; 37463c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 37563c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 37663c07aadSStefano Zampini if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported"); 37763c07aadSStefano Zampini } 37863c07aadSStefano Zampini ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 37963c07aadSStefano Zampini 38063c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 38163c07aadSStefano Zampini hdiag = hypre_ParCSRMatrixDiag(parcsr); 38263c07aadSStefano Zampini hoffd = hypre_ParCSRMatrixOffd(parcsr); 38363c07aadSStefano Zampini m = hypre_CSRMatrixNumRows(hdiag); 38463c07aadSStefano Zampini n = hypre_CSRMatrixNumCols(hdiag); 38563c07aadSStefano Zampini dnnz = hypre_CSRMatrixNumNonzeros(hdiag); 38663c07aadSStefano Zampini onnz = hypre_CSRMatrixNumNonzeros(hoffd); 38763c07aadSStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 38863c07aadSStefano Zampini ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr); 38963c07aadSStefano Zampini ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr); 39063c07aadSStefano Zampini ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr); 39163c07aadSStefano Zampini } else { 39263c07aadSStefano Zampini PetscInt nr; 39363c07aadSStefano Zampini PetscBool done; 39463c07aadSStefano Zampini if (size > 1) { 39563c07aadSStefano Zampini Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); 39663c07aadSStefano Zampini 39763c07aadSStefano Zampini ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); 39863c07aadSStefano 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); 39963c07aadSStefano 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); 40063c07aadSStefano Zampini ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr); 40163c07aadSStefano Zampini } else { 40263c07aadSStefano Zampini ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); 40363c07aadSStefano Zampini if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m); 40463c07aadSStefano 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); 40563c07aadSStefano Zampini ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr); 40663c07aadSStefano Zampini } 40763c07aadSStefano Zampini } 40863c07aadSStefano Zampini ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 40963c07aadSStefano Zampini ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr); 41063c07aadSStefano Zampini ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr); 41163c07aadSStefano Zampini iptr = djj; 41263c07aadSStefano Zampini aptr = da; 41363c07aadSStefano Zampini for (i=0; i<m; i++) { 41463c07aadSStefano Zampini PetscInt nc = dii[i+1]-dii[i]; 41563c07aadSStefano Zampini ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr); 41663c07aadSStefano Zampini iptr += nc; 41763c07aadSStefano Zampini aptr += nc; 41863c07aadSStefano Zampini } 41963c07aadSStefano Zampini if (size > 1) { 42063c07aadSStefano Zampini PetscInt *offdj,*coffd; 42163c07aadSStefano Zampini 42263c07aadSStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 42363c07aadSStefano Zampini ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr); 42463c07aadSStefano Zampini ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr); 42563c07aadSStefano Zampini ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr); 42663c07aadSStefano Zampini } else { 42763c07aadSStefano Zampini Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); 42863c07aadSStefano Zampini PetscInt nr,hr = hypre_CSRMatrixNumRows(hoffd); 42963c07aadSStefano Zampini PetscBool done; 43063c07aadSStefano Zampini 43163c07aadSStefano Zampini ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr); 43263c07aadSStefano 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); 43363c07aadSStefano 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); 43463c07aadSStefano Zampini ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr); 43563c07aadSStefano Zampini } 43663c07aadSStefano Zampini ierr = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 43763c07aadSStefano Zampini offdj = hypre_CSRMatrixJ(hoffd); 43863c07aadSStefano Zampini coffd = hypre_ParCSRMatrixColMapOffd(parcsr); 43963c07aadSStefano Zampini for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]]; 44063c07aadSStefano Zampini ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr); 44163c07aadSStefano Zampini iptr = ojj; 44263c07aadSStefano Zampini aptr = oa; 44363c07aadSStefano Zampini for (i=0; i<m; i++) { 44463c07aadSStefano Zampini PetscInt nc = oii[i+1]-oii[i]; 44563c07aadSStefano Zampini ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr); 44663c07aadSStefano Zampini iptr += nc; 44763c07aadSStefano Zampini aptr += nc; 44863c07aadSStefano Zampini } 44963c07aadSStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 45063c07aadSStefano Zampini Mat_MPIAIJ *b; 45163c07aadSStefano Zampini Mat_SeqAIJ *d,*o; 45263c07aadSStefano Zampini ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii, 45363c07aadSStefano Zampini djj,da,oii,ojj,oa,B);CHKERRQ(ierr); 45463c07aadSStefano Zampini /* hack MPIAIJ */ 45563c07aadSStefano Zampini b = (Mat_MPIAIJ*)((*B)->data); 45663c07aadSStefano Zampini d = (Mat_SeqAIJ*)b->A->data; 45763c07aadSStefano Zampini o = (Mat_SeqAIJ*)b->B->data; 45863c07aadSStefano Zampini d->free_a = PETSC_TRUE; 45963c07aadSStefano Zampini d->free_ij = PETSC_TRUE; 46063c07aadSStefano Zampini o->free_a = PETSC_TRUE; 46163c07aadSStefano Zampini o->free_ij = PETSC_TRUE; 46263c07aadSStefano Zampini } 46363c07aadSStefano Zampini } else if (reuse != MAT_REUSE_MATRIX) { 46463c07aadSStefano Zampini Mat_SeqAIJ* b; 46563c07aadSStefano Zampini ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr); 46663c07aadSStefano Zampini /* hack SeqAIJ */ 46763c07aadSStefano Zampini b = (Mat_SeqAIJ*)((*B)->data); 46863c07aadSStefano Zampini b->free_a = PETSC_TRUE; 46963c07aadSStefano Zampini b->free_ij = PETSC_TRUE; 47063c07aadSStefano Zampini } 47163c07aadSStefano Zampini if (reuse == MAT_INPLACE_MATRIX) { 47263c07aadSStefano Zampini ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 47363c07aadSStefano Zampini } 47463c07aadSStefano Zampini PetscFunctionReturn(0); 47563c07aadSStefano Zampini } 47663c07aadSStefano Zampini 47763c07aadSStefano Zampini #undef __FUNCT__ 47863c07aadSStefano Zampini #define __FUNCT__ "MatMultTranspose_HYPRE" 479ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y) 48063c07aadSStefano Zampini { 48163c07aadSStefano Zampini PetscErrorCode ierr; 48263c07aadSStefano Zampini 48363c07aadSStefano Zampini PetscFunctionBegin; 48463c07aadSStefano Zampini ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr); 48563c07aadSStefano Zampini PetscFunctionReturn(0); 48663c07aadSStefano Zampini } 48763c07aadSStefano Zampini 48863c07aadSStefano Zampini #undef __FUNCT__ 48963c07aadSStefano Zampini #define __FUNCT__ "MatMult_HYPRE" 490ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y) 49163c07aadSStefano Zampini { 49263c07aadSStefano Zampini PetscErrorCode ierr; 49363c07aadSStefano Zampini 49463c07aadSStefano Zampini PetscFunctionBegin; 49563c07aadSStefano Zampini ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr); 49663c07aadSStefano Zampini PetscFunctionReturn(0); 49763c07aadSStefano Zampini } 49863c07aadSStefano Zampini 49963c07aadSStefano Zampini #undef __FUNCT__ 50063c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_MultKernel_Private" 501ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans) 50263c07aadSStefano Zampini { 50363c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 50463c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 50563c07aadSStefano Zampini hypre_ParVector *hx,*hy; 50663c07aadSStefano Zampini PetscScalar *ax,*ay,*sax,*say; 50763c07aadSStefano Zampini PetscErrorCode ierr; 50863c07aadSStefano Zampini 50963c07aadSStefano Zampini PetscFunctionBegin; 51063c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 51163c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx)); 51263c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy)); 51363c07aadSStefano Zampini ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); 51463c07aadSStefano Zampini ierr = VecGetArray(y,&ay);CHKERRQ(ierr); 51563c07aadSStefano Zampini if (trans) { 51663c07aadSStefano Zampini HYPREReplacePointer(hA->x,ay,say); 51763c07aadSStefano Zampini HYPREReplacePointer(hA->b,ax,sax); 51863c07aadSStefano Zampini hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx); 51963c07aadSStefano Zampini HYPREReplacePointer(hA->x,say,ay); 52063c07aadSStefano Zampini HYPREReplacePointer(hA->b,sax,ax); 52163c07aadSStefano Zampini } else { 52263c07aadSStefano Zampini HYPREReplacePointer(hA->x,ax,sax); 52363c07aadSStefano Zampini HYPREReplacePointer(hA->b,ay,say); 52463c07aadSStefano Zampini hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy); 52563c07aadSStefano Zampini HYPREReplacePointer(hA->x,sax,ax); 52663c07aadSStefano Zampini HYPREReplacePointer(hA->b,say,ay); 52763c07aadSStefano Zampini } 52863c07aadSStefano Zampini ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); 52963c07aadSStefano Zampini ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr); 53063c07aadSStefano Zampini PetscFunctionReturn(0); 53163c07aadSStefano Zampini } 53263c07aadSStefano Zampini 53363c07aadSStefano Zampini #undef __FUNCT__ 53463c07aadSStefano Zampini #define __FUNCT__ "MatDestroy_HYPRE" 535ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A) 53663c07aadSStefano Zampini { 53763c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 53863c07aadSStefano Zampini PetscErrorCode ierr; 53963c07aadSStefano Zampini 54063c07aadSStefano Zampini PetscFunctionBegin; 54163c07aadSStefano Zampini if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x)); 54263c07aadSStefano Zampini if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b)); 543978814f1SStefano Zampini if (hA->ij) { 544978814f1SStefano Zampini if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL; 545978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij)); 546978814f1SStefano Zampini } 54763c07aadSStefano Zampini if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); } 54863c07aadSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr); 549*2df22349SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr); 55063c07aadSStefano Zampini ierr = PetscFree(A->data);CHKERRQ(ierr); 55163c07aadSStefano Zampini PetscFunctionReturn(0); 55263c07aadSStefano Zampini } 55363c07aadSStefano Zampini 55463c07aadSStefano Zampini #undef __FUNCT__ 55563c07aadSStefano Zampini #define __FUNCT__ "MatSetUp_HYPRE" 556ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A) 55763c07aadSStefano Zampini { 55863c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 55963c07aadSStefano Zampini Vec x,b; 56063c07aadSStefano Zampini PetscErrorCode ierr; 56163c07aadSStefano Zampini 56263c07aadSStefano Zampini PetscFunctionBegin; 56363c07aadSStefano Zampini if (hA->x) PetscFunctionReturn(0); 56463c07aadSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 56563c07aadSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 56663c07aadSStefano Zampini ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr); 56763c07aadSStefano Zampini ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr); 56863c07aadSStefano Zampini ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr); 56963c07aadSStefano Zampini ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr); 57063c07aadSStefano Zampini ierr = VecDestroy(&x);CHKERRQ(ierr); 57163c07aadSStefano Zampini ierr = VecDestroy(&b);CHKERRQ(ierr); 57263c07aadSStefano Zampini PetscFunctionReturn(0); 57363c07aadSStefano Zampini } 57463c07aadSStefano Zampini 57563c07aadSStefano Zampini #undef __FUNCT__ 57663c07aadSStefano Zampini #define __FUNCT__ "MatAssemblyEnd_HYPRE" 577ea9daf28SStefano Zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode) 57863c07aadSStefano Zampini { 57963c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 58063c07aadSStefano Zampini PetscErrorCode ierr; 58163c07aadSStefano Zampini 58263c07aadSStefano Zampini PetscFunctionBegin; 58363c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij)); 58463c07aadSStefano Zampini PetscFunctionReturn(0); 58563c07aadSStefano Zampini } 58663c07aadSStefano Zampini 587978814f1SStefano Zampini #undef __FUNCT__ 588978814f1SStefano Zampini #define __FUNCT__ "MatHYPRECreateFromParCSR" 589978814f1SStefano Zampini PETSC_EXTERN PetscErrorCode MatHYPRECreateFromParCSR(void *vparcsr, PetscCopyMode copymode, Mat* A) 590978814f1SStefano Zampini { 591978814f1SStefano Zampini Mat_HYPRE *hA; 59263c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 593978814f1SStefano Zampini MPI_Comm comm; 594978814f1SStefano Zampini PetscInt rstart,rend,cstart,cend,M,N; 595978814f1SStefano Zampini PetscErrorCode ierr; 596978814f1SStefano Zampini 597978814f1SStefano Zampini PetscFunctionBegin; 598978814f1SStefano Zampini parcsr = (hypre_ParCSRMatrix *)vparcsr; 599978814f1SStefano Zampini comm = hypre_ParCSRMatrixComm(parcsr); 600978814f1SStefano Zampini if (copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES"); 601978814f1SStefano Zampini 602978814f1SStefano Zampini /* access ParCSRMatrix */ 603978814f1SStefano Zampini rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr); 604978814f1SStefano Zampini rend = hypre_ParCSRMatrixLastRowIndex(parcsr); 605978814f1SStefano Zampini cstart = hypre_ParCSRMatrixFirstColDiag(parcsr); 606978814f1SStefano Zampini cend = hypre_ParCSRMatrixLastColDiag(parcsr); 607978814f1SStefano Zampini M = hypre_ParCSRMatrixGlobalNumRows(parcsr); 608978814f1SStefano Zampini N = hypre_ParCSRMatrixGlobalNumCols(parcsr); 609978814f1SStefano Zampini 610978814f1SStefano Zampini /* create PETSc matrix with MatHYPRE */ 611978814f1SStefano Zampini ierr = MatCreate(comm,A);CHKERRQ(ierr); 612978814f1SStefano Zampini ierr = MatSetSizes(*A,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr); 613978814f1SStefano Zampini ierr = MatSetType(*A,MATHYPRE);CHKERRQ(ierr); 614978814f1SStefano Zampini ierr = MatSetUp(*A);CHKERRQ(ierr); 615978814f1SStefano Zampini hA = (Mat_HYPRE*)((*A)->data); 616978814f1SStefano Zampini 617978814f1SStefano Zampini /* create HYPRE_IJMatrix */ 618978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij)); 619978814f1SStefano Zampini 620978814f1SStefano Zampini /* set ParCSR object */ 621978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); 622978814f1SStefano Zampini hypre_IJMatrixObject(hA->ij) = parcsr; 623978814f1SStefano Zampini 624978814f1SStefano Zampini /* set assembled flag */ 625978814f1SStefano Zampini hypre_IJMatrixAssembleFlag(hA->ij) = 1; 626978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij)); 627978814f1SStefano Zampini 628978814f1SStefano Zampini /* prevent from freeing the pointer */ 629978814f1SStefano Zampini if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE; 630978814f1SStefano Zampini 631978814f1SStefano Zampini ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 632978814f1SStefano Zampini ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 633978814f1SStefano Zampini PetscFunctionReturn(0); 634978814f1SStefano Zampini } 635978814f1SStefano Zampini 63663c07aadSStefano Zampini #undef __FUNCT__ 63763c07aadSStefano Zampini #define __FUNCT__ "MatCreate_HYPRE" 63863c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B) 63963c07aadSStefano Zampini { 64063c07aadSStefano Zampini Mat_HYPRE *hB; 64163c07aadSStefano Zampini PetscErrorCode ierr; 64263c07aadSStefano Zampini 64363c07aadSStefano Zampini PetscFunctionBegin; 64463c07aadSStefano Zampini ierr = PetscNewLog(B,&hB);CHKERRQ(ierr); 645978814f1SStefano Zampini hB->inner_free = PETSC_TRUE; 646978814f1SStefano Zampini 64763c07aadSStefano Zampini B->data = (void*)hB; 64863c07aadSStefano Zampini B->rmap->bs = 1; 64963c07aadSStefano Zampini B->assembled = PETSC_FALSE; 65063c07aadSStefano Zampini 65163c07aadSStefano Zampini B->ops->mult = MatMult_HYPRE; 65263c07aadSStefano Zampini B->ops->multtranspose = MatMultTranspose_HYPRE; 65363c07aadSStefano Zampini B->ops->setup = MatSetUp_HYPRE; 65463c07aadSStefano Zampini B->ops->destroy = MatDestroy_HYPRE; 65563c07aadSStefano Zampini B->ops->assemblyend = MatAssemblyEnd_HYPRE; 65663c07aadSStefano Zampini 65763c07aadSStefano Zampini ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr); 65863c07aadSStefano Zampini ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr); 65963c07aadSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr); 660*2df22349SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr); 66163c07aadSStefano Zampini PetscFunctionReturn(0); 66263c07aadSStefano Zampini } 66363c07aadSStefano Zampini 66463c07aadSStefano Zampini #if 0 66563c07aadSStefano Zampini /* 66663c07aadSStefano Zampini Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format 66763c07aadSStefano Zampini 66863c07aadSStefano Zampini This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first 66963c07aadSStefano Zampini which will corrupt the PETSc data structure if we did this. Need a work around to this problem. 67063c07aadSStefano Zampini */ 67163c07aadSStefano Zampini #include <_hypre_IJ_mv.h> 67263c07aadSStefano Zampini #include <HYPRE_IJ_mv.h> 67363c07aadSStefano Zampini 67463c07aadSStefano Zampini #undef __FUNCT__ 67563c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixLink" 67663c07aadSStefano Zampini PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij) 67763c07aadSStefano Zampini { 67863c07aadSStefano Zampini PetscErrorCode ierr; 67963c07aadSStefano Zampini PetscInt rstart,rend,cstart,cend; 68063c07aadSStefano Zampini PetscBool flg; 68163c07aadSStefano Zampini hypre_AuxParCSRMatrix *aux_matrix; 68263c07aadSStefano Zampini 68363c07aadSStefano Zampini PetscFunctionBegin; 68463c07aadSStefano Zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 68563c07aadSStefano Zampini PetscValidType(A,1); 68663c07aadSStefano Zampini PetscValidPointer(ij,2); 68763c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 68863c07aadSStefano Zampini if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices"); 68963c07aadSStefano Zampini ierr = MatSetUp(A);CHKERRQ(ierr); 69063c07aadSStefano Zampini 69163c07aadSStefano Zampini ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 69263c07aadSStefano Zampini rstart = A->rmap->rstart; 69363c07aadSStefano Zampini rend = A->rmap->rend; 69463c07aadSStefano Zampini cstart = A->cmap->rstart; 69563c07aadSStefano Zampini cend = A->cmap->rend; 69663c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij)); 69763c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR)); 69863c07aadSStefano Zampini 69963c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij)); 70063c07aadSStefano Zampini PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij)); 70163c07aadSStefano Zampini 70263c07aadSStefano Zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 70363c07aadSStefano Zampini 70463c07aadSStefano Zampini /* this is the Hack part where we monkey directly with the hypre datastructures */ 70563c07aadSStefano Zampini 70663c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij)); 70763c07aadSStefano Zampini ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr); 70863c07aadSStefano Zampini PetscFunctionReturn(0); 70963c07aadSStefano Zampini } 71063c07aadSStefano Zampini #endif 711