163c07aadSStefano Zampini 263c07aadSStefano Zampini /* 363c07aadSStefano Zampini Creates hypre ijmatrix from PETSc matrix 463c07aadSStefano Zampini */ 5225daaf8SStefano Zampini 6dd9c0a25Sstefano_zampini #include <petscmathypre.h> 763c07aadSStefano Zampini #include <petsc/private/matimpl.h> 863c07aadSStefano Zampini #include <../src/mat/impls/hypre/mhypre.h> 963c07aadSStefano Zampini #include <../src/mat/impls/aij/mpi/mpiaij.h> 1058968eb6SStefano Zampini #include <../src/vec/vec/impls/hypre/vhyp.h> 1158968eb6SStefano Zampini #include <HYPRE.h> 12c1a070e6SStefano Zampini #include <HYPRE_utilities.h> 13cd8bc7baSStefano Zampini #include <_hypre_parcsr_ls.h> 1468ec7858SStefano Zampini #include <_hypre_sstruct_ls.h> 1563c07aadSStefano Zampini 1675d48cdbSStefano Zampini PETSC_INTERN PetscErrorCode MatPtAP_IS_XAIJ(Mat,Mat,MatReuse,PetscReal,Mat*); 1775d48cdbSStefano Zampini 1863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*); 1963c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix); 2063c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix); 2163c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix); 22414bd5c3SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,PetscScalar,Vec,PetscScalar,Vec,PetscBool); 23225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void*); 24c69f721fSFande Kong PetscErrorCode MatSetValues_HYPRE(Mat, PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode ins); 2563c07aadSStefano Zampini 2663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij) 2763c07aadSStefano Zampini { 2863c07aadSStefano Zampini PetscErrorCode ierr; 2963c07aadSStefano Zampini PetscInt i,n_d,n_o; 3063c07aadSStefano Zampini const PetscInt *ia_d,*ia_o; 3163c07aadSStefano Zampini PetscBool done_d=PETSC_FALSE,done_o=PETSC_FALSE; 3263c07aadSStefano Zampini PetscInt *nnz_d=NULL,*nnz_o=NULL; 3363c07aadSStefano Zampini 3463c07aadSStefano Zampini PetscFunctionBegin; 3563c07aadSStefano Zampini if (A_d) { /* determine number of nonzero entries in local diagonal part */ 3663c07aadSStefano Zampini ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr); 3763c07aadSStefano Zampini if (done_d) { 3863c07aadSStefano Zampini ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr); 3963c07aadSStefano Zampini for (i=0; i<n_d; i++) { 4063c07aadSStefano Zampini nnz_d[i] = ia_d[i+1] - ia_d[i]; 4163c07aadSStefano Zampini } 4263c07aadSStefano Zampini } 4363c07aadSStefano Zampini ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr); 4463c07aadSStefano Zampini } 4563c07aadSStefano Zampini if (A_o) { /* determine number of nonzero entries in local off-diagonal part */ 4663c07aadSStefano Zampini ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr); 4763c07aadSStefano Zampini if (done_o) { 4863c07aadSStefano Zampini ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr); 4963c07aadSStefano Zampini for (i=0; i<n_o; i++) { 5063c07aadSStefano Zampini nnz_o[i] = ia_o[i+1] - ia_o[i]; 5163c07aadSStefano Zampini } 5263c07aadSStefano Zampini } 5363c07aadSStefano Zampini ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr); 5463c07aadSStefano Zampini } 5563c07aadSStefano Zampini if (done_d) { /* set number of nonzeros in HYPRE IJ matrix */ 5663c07aadSStefano Zampini if (!done_o) { /* only diagonal part */ 5763c07aadSStefano Zampini ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr); 5863c07aadSStefano Zampini for (i=0; i<n_d; i++) { 5963c07aadSStefano Zampini nnz_o[i] = 0; 6063c07aadSStefano Zampini } 6163c07aadSStefano Zampini } 6263c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o)); 6363c07aadSStefano Zampini ierr = PetscFree(nnz_d);CHKERRQ(ierr); 6463c07aadSStefano Zampini ierr = PetscFree(nnz_o);CHKERRQ(ierr); 6563c07aadSStefano Zampini } 6663c07aadSStefano Zampini PetscFunctionReturn(0); 6763c07aadSStefano Zampini } 6863c07aadSStefano Zampini 6963c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA) 7063c07aadSStefano Zampini { 7163c07aadSStefano Zampini PetscErrorCode ierr; 7263c07aadSStefano Zampini PetscInt rstart,rend,cstart,cend; 7363c07aadSStefano Zampini 7463c07aadSStefano Zampini PetscFunctionBegin; 75d975228cSstefano_zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 76d975228cSstefano_zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 7763c07aadSStefano Zampini rstart = A->rmap->rstart; 7863c07aadSStefano Zampini rend = A->rmap->rend; 7963c07aadSStefano Zampini cstart = A->cmap->rstart; 8063c07aadSStefano Zampini cend = A->cmap->rend; 8163c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij)); 8263c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); 8363c07aadSStefano Zampini { 8463c07aadSStefano Zampini PetscBool same; 8563c07aadSStefano Zampini Mat A_d,A_o; 8663c07aadSStefano Zampini const PetscInt *colmap; 8763c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr); 8863c07aadSStefano Zampini if (same) { 8963c07aadSStefano Zampini ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); 9063c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr); 9163c07aadSStefano Zampini PetscFunctionReturn(0); 9263c07aadSStefano Zampini } 9363c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr); 9463c07aadSStefano Zampini if (same) { 9563c07aadSStefano Zampini ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr); 9663c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr); 9763c07aadSStefano Zampini PetscFunctionReturn(0); 9863c07aadSStefano Zampini } 9963c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr); 10063c07aadSStefano Zampini if (same) { 10163c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr); 10263c07aadSStefano Zampini PetscFunctionReturn(0); 10363c07aadSStefano Zampini } 10463c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr); 10563c07aadSStefano Zampini if (same) { 10663c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr); 10763c07aadSStefano Zampini PetscFunctionReturn(0); 10863c07aadSStefano Zampini } 10963c07aadSStefano Zampini } 11063c07aadSStefano Zampini PetscFunctionReturn(0); 11163c07aadSStefano Zampini } 11263c07aadSStefano Zampini 11363c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij) 11463c07aadSStefano Zampini { 11563c07aadSStefano Zampini PetscErrorCode ierr; 11663c07aadSStefano Zampini PetscInt i,rstart,rend,ncols,nr,nc; 11763c07aadSStefano Zampini const PetscScalar *values; 11863c07aadSStefano Zampini const PetscInt *cols; 11963c07aadSStefano Zampini PetscBool flg; 12063c07aadSStefano Zampini 12163c07aadSStefano Zampini PetscFunctionBegin; 12263c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr); 12363c07aadSStefano Zampini ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr); 12463c07aadSStefano Zampini if (flg && nr == nc) { 12563c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr); 12663c07aadSStefano Zampini PetscFunctionReturn(0); 12763c07aadSStefano Zampini } 12863c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr); 12963c07aadSStefano Zampini if (flg) { 13063c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr); 13163c07aadSStefano Zampini PetscFunctionReturn(0); 13263c07aadSStefano Zampini } 13363c07aadSStefano Zampini 13463c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 13563c07aadSStefano Zampini ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 13663c07aadSStefano Zampini for (i=rstart; i<rend; i++) { 13763c07aadSStefano Zampini ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 138e3977e59Sstefano_zampini if (ncols) { 13963c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values)); 140e3977e59Sstefano_zampini } 14163c07aadSStefano Zampini ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr); 14263c07aadSStefano Zampini } 14363c07aadSStefano Zampini PetscFunctionReturn(0); 14463c07aadSStefano Zampini } 14563c07aadSStefano Zampini 14663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij) 14763c07aadSStefano Zampini { 14863c07aadSStefano Zampini PetscErrorCode ierr; 14963c07aadSStefano Zampini Mat_SeqAIJ *pdiag = (Mat_SeqAIJ*)A->data; 15058968eb6SStefano Zampini HYPRE_Int type; 15163c07aadSStefano Zampini hypre_ParCSRMatrix *par_matrix; 15263c07aadSStefano Zampini hypre_AuxParCSRMatrix *aux_matrix; 15363c07aadSStefano Zampini hypre_CSRMatrix *hdiag; 15463c07aadSStefano Zampini 15563c07aadSStefano Zampini PetscFunctionBegin; 15663c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij)); 157ea9daf28SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type)); 158ea9daf28SStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 159ea9daf28SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix)); 16063c07aadSStefano Zampini hdiag = hypre_ParCSRMatrixDiag(par_matrix); 16163c07aadSStefano Zampini /* 16263c07aadSStefano Zampini this is the Hack part where we monkey directly with the hypre datastructures 16363c07aadSStefano Zampini */ 1641d4906efSStefano Zampini ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));CHKERRQ(ierr); 1651d4906efSStefano Zampini ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));CHKERRQ(ierr); 1661d4906efSStefano Zampini ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));CHKERRQ(ierr); 167ea9daf28SStefano Zampini 168ea9daf28SStefano Zampini aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij); 16963c07aadSStefano Zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0; 17063c07aadSStefano Zampini PetscFunctionReturn(0); 17163c07aadSStefano Zampini } 17263c07aadSStefano Zampini 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 2192df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B) 2202df22349SStefano Zampini { 2212df22349SStefano Zampini Mat_HYPRE* mhA = (Mat_HYPRE*)(A->data); 2222df22349SStefano Zampini Mat lA; 2232df22349SStefano Zampini ISLocalToGlobalMapping rl2g,cl2g; 2242df22349SStefano Zampini IS is; 2252df22349SStefano Zampini hypre_ParCSRMatrix *hA; 2262df22349SStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 2272df22349SStefano Zampini MPI_Comm comm; 2282df22349SStefano Zampini PetscScalar *hdd,*hod,*aa,*data; 2297d968826Sstefano_zampini HYPRE_Int *col_map_offd,*hdi,*hdj,*hoi,*hoj; 2302df22349SStefano Zampini PetscInt *ii,*jj,*iptr,*jptr; 2312df22349SStefano Zampini PetscInt cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N; 23258968eb6SStefano Zampini HYPRE_Int type; 2332df22349SStefano Zampini PetscErrorCode ierr; 2342df22349SStefano Zampini 2352df22349SStefano Zampini PetscFunctionBegin; 236a1787963SStefano Zampini comm = PetscObjectComm((PetscObject)A); 2372df22349SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type)); 2382df22349SStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 2392df22349SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA)); 2402df22349SStefano Zampini M = hypre_ParCSRMatrixGlobalNumRows(hA); 2412df22349SStefano Zampini N = hypre_ParCSRMatrixGlobalNumCols(hA); 2422df22349SStefano Zampini str = hypre_ParCSRMatrixFirstRowIndex(hA); 2432df22349SStefano Zampini stc = hypre_ParCSRMatrixFirstColDiag(hA); 2442df22349SStefano Zampini hdiag = hypre_ParCSRMatrixDiag(hA); 2452df22349SStefano Zampini hoffd = hypre_ParCSRMatrixOffd(hA); 2462df22349SStefano Zampini dr = hypre_CSRMatrixNumRows(hdiag); 2472df22349SStefano Zampini dc = hypre_CSRMatrixNumCols(hdiag); 2482df22349SStefano Zampini nnz = hypre_CSRMatrixNumNonzeros(hdiag); 2492df22349SStefano Zampini hdi = hypre_CSRMatrixI(hdiag); 2502df22349SStefano Zampini hdj = hypre_CSRMatrixJ(hdiag); 2512df22349SStefano Zampini hdd = hypre_CSRMatrixData(hdiag); 2522df22349SStefano Zampini oc = hypre_CSRMatrixNumCols(hoffd); 2532df22349SStefano Zampini nnz += hypre_CSRMatrixNumNonzeros(hoffd); 2542df22349SStefano Zampini hoi = hypre_CSRMatrixI(hoffd); 2552df22349SStefano Zampini hoj = hypre_CSRMatrixJ(hoffd); 2562df22349SStefano Zampini hod = hypre_CSRMatrixData(hoffd); 2572df22349SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 2582df22349SStefano Zampini PetscInt *aux; 2592df22349SStefano Zampini 2602df22349SStefano Zampini /* generate l2g maps for rows and cols */ 2612df22349SStefano Zampini ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr); 2622df22349SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr); 2632df22349SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 2642df22349SStefano Zampini col_map_offd = hypre_ParCSRMatrixColMapOffd(hA); 2652df22349SStefano Zampini ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr); 2662df22349SStefano Zampini for (i=0; i<dc; i++) aux[i] = i+stc; 2672df22349SStefano Zampini for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i]; 2682df22349SStefano Zampini ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr); 2692df22349SStefano Zampini ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr); 2702df22349SStefano Zampini ierr = ISDestroy(&is);CHKERRQ(ierr); 2712df22349SStefano Zampini /* create MATIS object */ 2722df22349SStefano Zampini ierr = MatCreate(comm,B);CHKERRQ(ierr); 2732df22349SStefano Zampini ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr); 2742df22349SStefano Zampini ierr = MatSetType(*B,MATIS);CHKERRQ(ierr); 2752df22349SStefano Zampini ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr); 2762df22349SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr); 2772df22349SStefano Zampini ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr); 2782df22349SStefano Zampini 2792df22349SStefano Zampini /* allocate CSR for local matrix */ 2802df22349SStefano Zampini ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr); 2812df22349SStefano Zampini ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr); 2822df22349SStefano Zampini ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr); 2832df22349SStefano Zampini } else { 2842df22349SStefano Zampini PetscInt nr; 2852df22349SStefano Zampini PetscBool done; 2862df22349SStefano Zampini ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr); 2872df22349SStefano Zampini ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr); 2882df22349SStefano 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); 2892df22349SStefano 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); 2902df22349SStefano Zampini ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr); 2912df22349SStefano Zampini } 2922df22349SStefano Zampini /* merge local matrices */ 2932df22349SStefano Zampini ii = iptr; 2942df22349SStefano Zampini jj = jptr; 2952df22349SStefano Zampini aa = data; 2962df22349SStefano Zampini *ii = *(hdi++) + *(hoi++); 2972df22349SStefano Zampini for (jd=0,jo=0,cum=0; *ii<nnz; cum++) { 2982df22349SStefano Zampini PetscScalar *aold = aa; 2992df22349SStefano Zampini PetscInt *jold = jj,nc = jd+jo; 3002df22349SStefano Zampini for (; jd<*hdi; jd++) { *jj++ = *hdj++; *aa++ = *hdd++; } 3012df22349SStefano Zampini for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; } 3022df22349SStefano Zampini *(++ii) = *(hdi++) + *(hoi++); 3032df22349SStefano Zampini ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr); 3042df22349SStefano Zampini } 3052df22349SStefano Zampini for (; cum<dr; cum++) *(++ii) = nnz; 3062df22349SStefano Zampini if (reuse != MAT_REUSE_MATRIX) { 307a033916dSStefano Zampini Mat_SeqAIJ* a; 308a033916dSStefano Zampini 3092df22349SStefano Zampini ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr); 3102df22349SStefano Zampini ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr); 311a033916dSStefano Zampini /* hack SeqAIJ */ 312a033916dSStefano Zampini a = (Mat_SeqAIJ*)(lA->data); 313a033916dSStefano Zampini a->free_a = PETSC_TRUE; 314a033916dSStefano Zampini a->free_ij = PETSC_TRUE; 3152df22349SStefano Zampini ierr = MatDestroy(&lA);CHKERRQ(ierr); 3162df22349SStefano Zampini } 3172df22349SStefano Zampini ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3182df22349SStefano Zampini ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 3192df22349SStefano Zampini if (reuse == MAT_INPLACE_MATRIX) { 3202df22349SStefano Zampini ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 3212df22349SStefano Zampini } 3222df22349SStefano Zampini PetscFunctionReturn(0); 3232df22349SStefano Zampini } 3242df22349SStefano Zampini 32563c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B) 32663c07aadSStefano Zampini { 32784d4e069SStefano Zampini Mat M = NULL; 32863c07aadSStefano Zampini Mat_HYPRE *hB; 32963c07aadSStefano Zampini MPI_Comm comm = PetscObjectComm((PetscObject)A); 33063c07aadSStefano Zampini PetscErrorCode ierr; 33163c07aadSStefano Zampini 33263c07aadSStefano Zampini PetscFunctionBegin; 33363c07aadSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 33463c07aadSStefano Zampini /* always destroy the old matrix and create a new memory; 33563c07aadSStefano Zampini hope this does not churn the memory too much. The problem 33663c07aadSStefano Zampini is I do not know if it is possible to put the matrix back to 33763c07aadSStefano Zampini its initial state so that we can directly copy the values 33863c07aadSStefano Zampini the second time through. */ 33963c07aadSStefano Zampini hB = (Mat_HYPRE*)((*B)->data); 34063c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij)); 34163c07aadSStefano Zampini } else { 34284d4e069SStefano Zampini ierr = MatCreate(comm,&M);CHKERRQ(ierr); 34384d4e069SStefano Zampini ierr = MatSetType(M,MATHYPRE);CHKERRQ(ierr); 34484d4e069SStefano Zampini ierr = MatSetSizes(M,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr); 34584d4e069SStefano Zampini hB = (Mat_HYPRE*)(M->data); 34684d4e069SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) *B = M; 34763c07aadSStefano Zampini } 34863c07aadSStefano Zampini ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr); 34963c07aadSStefano Zampini ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr); 35084d4e069SStefano Zampini if (reuse == MAT_INPLACE_MATRIX) { 35184d4e069SStefano Zampini ierr = MatHeaderReplace(A,&M);CHKERRQ(ierr); 35284d4e069SStefano Zampini } 3534ec6421dSstefano_zampini (*B)->preallocated = PETSC_TRUE; 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 359ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B) 36063c07aadSStefano Zampini { 36163c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 36263c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 36363c07aadSStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 36463c07aadSStefano Zampini MPI_Comm comm; 36563c07aadSStefano Zampini PetscScalar *da,*oa,*aptr; 36663c07aadSStefano Zampini PetscInt *dii,*djj,*oii,*ojj,*iptr; 36763c07aadSStefano Zampini PetscInt i,dnnz,onnz,m,n; 36858968eb6SStefano Zampini HYPRE_Int type; 36963c07aadSStefano Zampini PetscMPIInt size; 37063c07aadSStefano Zampini PetscErrorCode ierr; 37163c07aadSStefano Zampini 37263c07aadSStefano Zampini PetscFunctionBegin; 37363c07aadSStefano Zampini comm = PetscObjectComm((PetscObject)A); 37463c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 37563c07aadSStefano Zampini if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 37663c07aadSStefano Zampini if (reuse == MAT_REUSE_MATRIX) { 37763c07aadSStefano Zampini PetscBool ismpiaij,isseqaij; 37863c07aadSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 3794099cc6bSBarry Smith ierr = PetscObjectBaseTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 38063c07aadSStefano Zampini if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported"); 38163c07aadSStefano Zampini } 38263c07aadSStefano Zampini ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 38363c07aadSStefano Zampini 38463c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 38563c07aadSStefano Zampini hdiag = hypre_ParCSRMatrixDiag(parcsr); 38663c07aadSStefano Zampini hoffd = hypre_ParCSRMatrixOffd(parcsr); 38763c07aadSStefano Zampini m = hypre_CSRMatrixNumRows(hdiag); 38863c07aadSStefano Zampini n = hypre_CSRMatrixNumCols(hdiag); 38963c07aadSStefano Zampini dnnz = hypre_CSRMatrixNumNonzeros(hdiag); 39063c07aadSStefano Zampini onnz = hypre_CSRMatrixNumNonzeros(hoffd); 391225daaf8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 39263c07aadSStefano Zampini ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr); 39363c07aadSStefano Zampini ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr); 39463c07aadSStefano Zampini ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr); 395225daaf8SStefano Zampini } else if (reuse == MAT_REUSE_MATRIX) { 39663c07aadSStefano Zampini PetscInt nr; 39763c07aadSStefano Zampini PetscBool done; 39863c07aadSStefano Zampini if (size > 1) { 39963c07aadSStefano Zampini Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); 40063c07aadSStefano Zampini 40163c07aadSStefano Zampini ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); 40263c07aadSStefano 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); 40363c07aadSStefano 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); 40463c07aadSStefano Zampini ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr); 40563c07aadSStefano Zampini } else { 40663c07aadSStefano Zampini ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr); 40763c07aadSStefano Zampini if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m); 40863c07aadSStefano 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); 40963c07aadSStefano Zampini ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr); 41063c07aadSStefano Zampini } 411225daaf8SStefano Zampini } else { /* MAT_INPLACE_MATRIX */ 4127d968826Sstefano_zampini dii = (PetscInt*)hypre_CSRMatrixI(hdiag); 4137d968826Sstefano_zampini djj = (PetscInt*)hypre_CSRMatrixJ(hdiag); 414225daaf8SStefano Zampini da = hypre_CSRMatrixData(hdiag); 41563c07aadSStefano Zampini } 41663c07aadSStefano Zampini ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 41763c07aadSStefano Zampini ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr); 41863c07aadSStefano Zampini ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr); 41963c07aadSStefano Zampini iptr = djj; 42063c07aadSStefano Zampini aptr = da; 42163c07aadSStefano Zampini for (i=0; i<m; i++) { 42263c07aadSStefano Zampini PetscInt nc = dii[i+1]-dii[i]; 42363c07aadSStefano Zampini ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr); 42463c07aadSStefano Zampini iptr += nc; 42563c07aadSStefano Zampini aptr += nc; 42663c07aadSStefano Zampini } 42763c07aadSStefano Zampini if (size > 1) { 4287d968826Sstefano_zampini HYPRE_Int *offdj,*coffd; 42963c07aadSStefano Zampini 430225daaf8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 43163c07aadSStefano Zampini ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr); 43263c07aadSStefano Zampini ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr); 43363c07aadSStefano Zampini ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr); 434225daaf8SStefano Zampini } else if (reuse == MAT_REUSE_MATRIX) { 43563c07aadSStefano Zampini Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data); 43663c07aadSStefano Zampini PetscInt nr,hr = hypre_CSRMatrixNumRows(hoffd); 43763c07aadSStefano Zampini PetscBool done; 43863c07aadSStefano Zampini 43963c07aadSStefano Zampini ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr); 44063c07aadSStefano 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); 44163c07aadSStefano 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); 44263c07aadSStefano Zampini ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr); 443225daaf8SStefano Zampini } else { /* MAT_INPLACE_MATRIX */ 4447d968826Sstefano_zampini oii = (PetscInt*)hypre_CSRMatrixI(hoffd); 4457d968826Sstefano_zampini ojj = (PetscInt*)hypre_CSRMatrixJ(hoffd); 446225daaf8SStefano Zampini oa = hypre_CSRMatrixData(hoffd); 44763c07aadSStefano Zampini } 44863c07aadSStefano Zampini ierr = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr); 44963c07aadSStefano Zampini offdj = hypre_CSRMatrixJ(hoffd); 45063c07aadSStefano Zampini coffd = hypre_ParCSRMatrixColMapOffd(parcsr); 45163c07aadSStefano Zampini for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]]; 45263c07aadSStefano Zampini ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr); 45363c07aadSStefano Zampini iptr = ojj; 45463c07aadSStefano Zampini aptr = oa; 45563c07aadSStefano Zampini for (i=0; i<m; i++) { 45663c07aadSStefano Zampini PetscInt nc = oii[i+1]-oii[i]; 45763c07aadSStefano Zampini ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr); 45863c07aadSStefano Zampini iptr += nc; 45963c07aadSStefano Zampini aptr += nc; 46063c07aadSStefano Zampini } 461225daaf8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 46263c07aadSStefano Zampini Mat_MPIAIJ *b; 46363c07aadSStefano Zampini Mat_SeqAIJ *d,*o; 464225daaf8SStefano Zampini 46541073d71Sstefano_zampini ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,B);CHKERRQ(ierr); 46663c07aadSStefano Zampini /* hack MPIAIJ */ 46763c07aadSStefano Zampini b = (Mat_MPIAIJ*)((*B)->data); 46863c07aadSStefano Zampini d = (Mat_SeqAIJ*)b->A->data; 46963c07aadSStefano Zampini o = (Mat_SeqAIJ*)b->B->data; 47063c07aadSStefano Zampini d->free_a = PETSC_TRUE; 47163c07aadSStefano Zampini d->free_ij = PETSC_TRUE; 47263c07aadSStefano Zampini o->free_a = PETSC_TRUE; 47363c07aadSStefano Zampini o->free_ij = PETSC_TRUE; 474225daaf8SStefano Zampini } else if (reuse == MAT_INPLACE_MATRIX) { 475225daaf8SStefano Zampini Mat T; 47641073d71Sstefano_zampini ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,&T);CHKERRQ(ierr); 477225daaf8SStefano Zampini hypre_CSRMatrixI(hdiag) = NULL; 478225daaf8SStefano Zampini hypre_CSRMatrixJ(hdiag) = NULL; 479225daaf8SStefano Zampini hypre_CSRMatrixData(hdiag) = NULL; 480225daaf8SStefano Zampini hypre_CSRMatrixI(hoffd) = NULL; 481225daaf8SStefano Zampini hypre_CSRMatrixJ(hoffd) = NULL; 482225daaf8SStefano Zampini hypre_CSRMatrixData(hoffd) = NULL; 483225daaf8SStefano Zampini ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr); 48463c07aadSStefano Zampini } 485225daaf8SStefano Zampini } else { 486225daaf8SStefano Zampini oii = NULL; 487225daaf8SStefano Zampini ojj = NULL; 488225daaf8SStefano Zampini oa = NULL; 489225daaf8SStefano Zampini if (reuse == MAT_INITIAL_MATRIX) { 49063c07aadSStefano Zampini Mat_SeqAIJ* b; 49163c07aadSStefano Zampini ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr); 49263c07aadSStefano Zampini /* hack SeqAIJ */ 49363c07aadSStefano Zampini b = (Mat_SeqAIJ*)((*B)->data); 49463c07aadSStefano Zampini b->free_a = PETSC_TRUE; 49563c07aadSStefano Zampini b->free_ij = PETSC_TRUE; 496225daaf8SStefano Zampini } else if (reuse == MAT_INPLACE_MATRIX) { 497225daaf8SStefano Zampini Mat T; 498225daaf8SStefano Zampini ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,&T);CHKERRQ(ierr); 499225daaf8SStefano Zampini hypre_CSRMatrixI(hdiag) = NULL; 500225daaf8SStefano Zampini hypre_CSRMatrixJ(hdiag) = NULL; 501225daaf8SStefano Zampini hypre_CSRMatrixData(hdiag) = NULL; 502225daaf8SStefano Zampini ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr); 50363c07aadSStefano Zampini } 504225daaf8SStefano Zampini } 505225daaf8SStefano Zampini 506225daaf8SStefano Zampini /* we have to use hypre_Tfree to free the arrays */ 50763c07aadSStefano Zampini if (reuse == MAT_INPLACE_MATRIX) { 508225daaf8SStefano Zampini void *ptrs[6] = {dii,djj,da,oii,ojj,oa}; 509225daaf8SStefano Zampini const char *names[6] = {"_hypre_csr_dii", 510225daaf8SStefano Zampini "_hypre_csr_djj", 511225daaf8SStefano Zampini "_hypre_csr_da", 512225daaf8SStefano Zampini "_hypre_csr_oii", 513225daaf8SStefano Zampini "_hypre_csr_ojj", 514225daaf8SStefano Zampini "_hypre_csr_oa"}; 515225daaf8SStefano Zampini for (i=0; i<6; i++) { 516225daaf8SStefano Zampini PetscContainer c; 517225daaf8SStefano Zampini 518225daaf8SStefano Zampini ierr = PetscContainerCreate(comm,&c);CHKERRQ(ierr); 519225daaf8SStefano Zampini ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr); 520225daaf8SStefano Zampini ierr = PetscContainerSetUserDestroy(c,hypre_array_destroy);CHKERRQ(ierr); 521225daaf8SStefano Zampini ierr = PetscObjectCompose((PetscObject)(*B),names[i],(PetscObject)c);CHKERRQ(ierr); 522225daaf8SStefano Zampini ierr = PetscContainerDestroy(&c);CHKERRQ(ierr); 523225daaf8SStefano Zampini } 52463c07aadSStefano Zampini } 52563c07aadSStefano Zampini PetscFunctionReturn(0); 52663c07aadSStefano Zampini } 52763c07aadSStefano Zampini 528613e5ff0Sstefano_zampini static PetscErrorCode MatAIJGetParCSR_Private(Mat A, hypre_ParCSRMatrix **hA) 529c1a070e6SStefano Zampini { 530613e5ff0Sstefano_zampini hypre_ParCSRMatrix *tA; 531c1a070e6SStefano Zampini hypre_CSRMatrix *hdiag,*hoffd; 532c1a070e6SStefano Zampini Mat_SeqAIJ *diag,*offd; 533c1a070e6SStefano Zampini PetscInt *garray,noffd,dnnz,onnz,*row_starts,*col_starts; 534c1a070e6SStefano Zampini MPI_Comm comm = PetscObjectComm((PetscObject)A); 535613e5ff0Sstefano_zampini PetscBool ismpiaij,isseqaij; 536c1a070e6SStefano Zampini PetscErrorCode ierr; 537c1a070e6SStefano Zampini 538c1a070e6SStefano Zampini PetscFunctionBegin; 539c1a070e6SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 5404099cc6bSBarry Smith ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 541c1a070e6SStefano Zampini if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type); 542c1a070e6SStefano Zampini if (ismpiaij) { 543c1a070e6SStefano Zampini Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data); 544c1a070e6SStefano Zampini 545c1a070e6SStefano Zampini diag = (Mat_SeqAIJ*)a->A->data; 546c1a070e6SStefano Zampini offd = (Mat_SeqAIJ*)a->B->data; 547c1a070e6SStefano Zampini garray = a->garray; 548c1a070e6SStefano Zampini noffd = a->B->cmap->N; 549c1a070e6SStefano Zampini dnnz = diag->nz; 550c1a070e6SStefano Zampini onnz = offd->nz; 551c1a070e6SStefano Zampini } else { 552c1a070e6SStefano Zampini diag = (Mat_SeqAIJ*)A->data; 553c1a070e6SStefano Zampini offd = NULL; 554c1a070e6SStefano Zampini garray = NULL; 555c1a070e6SStefano Zampini noffd = 0; 556c1a070e6SStefano Zampini dnnz = diag->nz; 557c1a070e6SStefano Zampini onnz = 0; 558c1a070e6SStefano Zampini } 559225daaf8SStefano Zampini 560c1a070e6SStefano Zampini /* create a temporary ParCSR */ 561c1a070e6SStefano Zampini if (HYPRE_AssumedPartitionCheck()) { 562c1a070e6SStefano Zampini PetscMPIInt myid; 563c1a070e6SStefano Zampini 564c1a070e6SStefano Zampini ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr); 565c1a070e6SStefano Zampini row_starts = A->rmap->range + myid; 566c1a070e6SStefano Zampini col_starts = A->cmap->range + myid; 567c1a070e6SStefano Zampini } else { 568c1a070e6SStefano Zampini row_starts = A->rmap->range; 569c1a070e6SStefano Zampini col_starts = A->cmap->range; 570c1a070e6SStefano Zampini } 5717d968826Sstefano_zampini tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,(HYPRE_Int*)row_starts,(HYPRE_Int*)col_starts,noffd,dnnz,onnz); 572c1a070e6SStefano Zampini hypre_ParCSRMatrixSetRowStartsOwner(tA,0); 573c1a070e6SStefano Zampini hypre_ParCSRMatrixSetColStartsOwner(tA,0); 574c1a070e6SStefano Zampini 575225daaf8SStefano Zampini /* set diagonal part */ 576c1a070e6SStefano Zampini hdiag = hypre_ParCSRMatrixDiag(tA); 5777d968826Sstefano_zampini hypre_CSRMatrixI(hdiag) = (HYPRE_Int*)diag->i; 5787d968826Sstefano_zampini hypre_CSRMatrixJ(hdiag) = (HYPRE_Int*)diag->j; 579c1a070e6SStefano Zampini hypre_CSRMatrixData(hdiag) = diag->a; 580c1a070e6SStefano Zampini hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz; 581c1a070e6SStefano Zampini hypre_CSRMatrixSetRownnz(hdiag); 582c1a070e6SStefano Zampini hypre_CSRMatrixSetDataOwner(hdiag,0); 583c1a070e6SStefano Zampini 584225daaf8SStefano Zampini /* set offdiagonal part */ 585c1a070e6SStefano Zampini hoffd = hypre_ParCSRMatrixOffd(tA); 586c1a070e6SStefano Zampini if (offd) { 5877d968826Sstefano_zampini hypre_CSRMatrixI(hoffd) = (HYPRE_Int*)offd->i; 5887d968826Sstefano_zampini hypre_CSRMatrixJ(hoffd) = (HYPRE_Int*)offd->j; 589c1a070e6SStefano Zampini hypre_CSRMatrixData(hoffd) = offd->a; 590c1a070e6SStefano Zampini hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz; 591c1a070e6SStefano Zampini hypre_CSRMatrixSetRownnz(hoffd); 592c1a070e6SStefano Zampini hypre_CSRMatrixSetDataOwner(hoffd,0); 593c1a070e6SStefano Zampini hypre_ParCSRMatrixSetNumNonzeros(tA); 5947d968826Sstefano_zampini hypre_ParCSRMatrixColMapOffd(tA) = (HYPRE_Int*)garray; 595c1a070e6SStefano Zampini } 596613e5ff0Sstefano_zampini *hA = tA; 597613e5ff0Sstefano_zampini PetscFunctionReturn(0); 598613e5ff0Sstefano_zampini } 599c1a070e6SStefano Zampini 600613e5ff0Sstefano_zampini static PetscErrorCode MatAIJRestoreParCSR_Private(Mat A, hypre_ParCSRMatrix **hA) 601613e5ff0Sstefano_zampini { 602613e5ff0Sstefano_zampini hypre_CSRMatrix *hdiag,*hoffd; 603c1a070e6SStefano Zampini 604613e5ff0Sstefano_zampini PetscFunctionBegin; 605613e5ff0Sstefano_zampini hdiag = hypre_ParCSRMatrixDiag(*hA); 606613e5ff0Sstefano_zampini hoffd = hypre_ParCSRMatrixOffd(*hA); 607c1a070e6SStefano Zampini /* set pointers to NULL before destroying tA */ 608c1a070e6SStefano Zampini hypre_CSRMatrixI(hdiag) = NULL; 609c1a070e6SStefano Zampini hypre_CSRMatrixJ(hdiag) = NULL; 610c1a070e6SStefano Zampini hypre_CSRMatrixData(hdiag) = NULL; 611c1a070e6SStefano Zampini hypre_CSRMatrixI(hoffd) = NULL; 612c1a070e6SStefano Zampini hypre_CSRMatrixJ(hoffd) = NULL; 613c1a070e6SStefano Zampini hypre_CSRMatrixData(hoffd) = NULL; 614613e5ff0Sstefano_zampini hypre_ParCSRMatrixColMapOffd(*hA) = NULL; 615613e5ff0Sstefano_zampini hypre_ParCSRMatrixDestroy(*hA); 616613e5ff0Sstefano_zampini *hA = NULL; 617613e5ff0Sstefano_zampini PetscFunctionReturn(0); 618613e5ff0Sstefano_zampini } 619613e5ff0Sstefano_zampini 620613e5ff0Sstefano_zampini /* calls RAP from BoomerAMG: 6213dad0653Sstefano_zampini the resulting ParCSR will not own the column and row starts 6223dad0653Sstefano_zampini It looks like we don't need to have the diagonal entries 6233dad0653Sstefano_zampini ordered first in the rows of the diagonal part 6243dad0653Sstefano_zampini for boomerAMGBuildCoarseOperator to work */ 625a055b5aaSBarry Smith static PetscErrorCode MatHYPRE_ParCSR_RAP(hypre_ParCSRMatrix *hR, hypre_ParCSRMatrix *hA,hypre_ParCSRMatrix *hP, hypre_ParCSRMatrix **hRAP) 626613e5ff0Sstefano_zampini { 627613e5ff0Sstefano_zampini HYPRE_Int P_owns_col_starts,R_owns_row_starts; 628613e5ff0Sstefano_zampini 629613e5ff0Sstefano_zampini PetscFunctionBegin; 630613e5ff0Sstefano_zampini P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(hP); 631613e5ff0Sstefano_zampini R_owns_row_starts = hypre_ParCSRMatrixOwnsRowStarts(hR); 632613e5ff0Sstefano_zampini PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(hR,hA,hP,hRAP)); 633613e5ff0Sstefano_zampini PetscStackCallStandard(hypre_ParCSRMatrixSetNumNonzeros,(*hRAP)); 634613e5ff0Sstefano_zampini /* hypre_BoomerAMGBuildCoarseOperator steals the col_starts from P and the row_starts from R */ 635613e5ff0Sstefano_zampini hypre_ParCSRMatrixSetRowStartsOwner(*hRAP,0); 636613e5ff0Sstefano_zampini hypre_ParCSRMatrixSetColStartsOwner(*hRAP,0); 637613e5ff0Sstefano_zampini if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(hP,1); 638613e5ff0Sstefano_zampini if (R_owns_row_starts) hypre_ParCSRMatrixSetRowStartsOwner(hR,1); 639613e5ff0Sstefano_zampini PetscFunctionReturn(0); 640613e5ff0Sstefano_zampini } 641613e5ff0Sstefano_zampini 6426f231fbdSstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat P,Mat C) 643613e5ff0Sstefano_zampini { 6446f231fbdSstefano_zampini Mat B; 645613e5ff0Sstefano_zampini hypre_ParCSRMatrix *hA,*hP,*hPtAP; 646613e5ff0Sstefano_zampini PetscErrorCode ierr; 647613e5ff0Sstefano_zampini 648613e5ff0Sstefano_zampini PetscFunctionBegin; 649613e5ff0Sstefano_zampini ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr); 650613e5ff0Sstefano_zampini ierr = MatAIJGetParCSR_Private(P,&hP);CHKERRQ(ierr); 651613e5ff0Sstefano_zampini ierr = MatHYPRE_ParCSR_RAP(hP,hA,hP,&hPtAP);CHKERRQ(ierr); 6526f231fbdSstefano_zampini ierr = MatCreateFromParCSR(hPtAP,MATAIJ,PETSC_OWN_POINTER,&B);CHKERRQ(ierr); 6536f231fbdSstefano_zampini ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr); 654613e5ff0Sstefano_zampini ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr); 655613e5ff0Sstefano_zampini ierr = MatAIJRestoreParCSR_Private(P,&hP);CHKERRQ(ierr); 6566f231fbdSstefano_zampini PetscFunctionReturn(0); 6576f231fbdSstefano_zampini } 6586f231fbdSstefano_zampini 6596f231fbdSstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat P,PetscReal fill,Mat *C) 6606f231fbdSstefano_zampini { 6616f231fbdSstefano_zampini PetscErrorCode ierr; 662ab4d48faSStefano Zampini 6636f231fbdSstefano_zampini PetscFunctionBegin; 6646f231fbdSstefano_zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr); 6656f231fbdSstefano_zampini ierr = MatSetType(*C,MATAIJ);CHKERRQ(ierr); 6666f231fbdSstefano_zampini (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_AIJ_wHYPRE; 667613e5ff0Sstefano_zampini PetscFunctionReturn(0); 668613e5ff0Sstefano_zampini } 669613e5ff0Sstefano_zampini 6704cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_HYPRE(Mat A,Mat P,Mat C) 671613e5ff0Sstefano_zampini { 6724cc28894Sstefano_zampini Mat B; 6734cc28894Sstefano_zampini Mat_HYPRE *hP; 674613e5ff0Sstefano_zampini hypre_ParCSRMatrix *hA = NULL,*Pparcsr,*ptapparcsr; 675613e5ff0Sstefano_zampini HYPRE_Int type; 676613e5ff0Sstefano_zampini MPI_Comm comm = PetscObjectComm((PetscObject)A); 6774cc28894Sstefano_zampini PetscBool ishypre; 678613e5ff0Sstefano_zampini PetscErrorCode ierr; 679613e5ff0Sstefano_zampini 680613e5ff0Sstefano_zampini PetscFunctionBegin; 6814cc28894Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr); 6824cc28894Sstefano_zampini if (!ishypre) SETERRQ1(comm,PETSC_ERR_USER,"P should be of type %s",MATHYPRE); 6834cc28894Sstefano_zampini hP = (Mat_HYPRE*)P->data; 684613e5ff0Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type)); 685613e5ff0Sstefano_zampini if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 686613e5ff0Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr)); 687613e5ff0Sstefano_zampini 688613e5ff0Sstefano_zampini ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr); 689613e5ff0Sstefano_zampini ierr = MatHYPRE_ParCSR_RAP(Pparcsr,hA,Pparcsr,&ptapparcsr);CHKERRQ(ierr); 690613e5ff0Sstefano_zampini ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr); 691225daaf8SStefano Zampini 6924cc28894Sstefano_zampini /* create temporary matrix and merge to C */ 6934cc28894Sstefano_zampini ierr = MatCreateFromParCSR(ptapparcsr,((PetscObject)C)->type_name,PETSC_OWN_POINTER,&B);CHKERRQ(ierr); 6944cc28894Sstefano_zampini ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr); 6954cc28894Sstefano_zampini PetscFunctionReturn(0); 6964cc28894Sstefano_zampini } 6974cc28894Sstefano_zampini 6984cc28894Sstefano_zampini static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 6994cc28894Sstefano_zampini { 7004cc28894Sstefano_zampini PetscErrorCode ierr; 7014cc28894Sstefano_zampini 7024cc28894Sstefano_zampini PetscFunctionBegin; 7034cc28894Sstefano_zampini if (scall == MAT_INITIAL_MATRIX) { 7044cc28894Sstefano_zampini const char *deft = MATAIJ; 7054cc28894Sstefano_zampini char type[256]; 7064cc28894Sstefano_zampini PetscBool flg; 7074cc28894Sstefano_zampini 7084cc28894Sstefano_zampini ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 7094cc28894Sstefano_zampini ierr = PetscOptionsFList("-matptap_hypre_outtype","Matrix type",NULL,MatList,deft,type,256,&flg);CHKERRQ(ierr); 7104cc28894Sstefano_zampini ierr = PetscOptionsEnd();CHKERRQ(ierr); 7114cc28894Sstefano_zampini ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 7124cc28894Sstefano_zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr); 7134cc28894Sstefano_zampini if (flg) { 7144cc28894Sstefano_zampini ierr = MatSetType(*C,type);CHKERRQ(ierr); 7154cc28894Sstefano_zampini } else { 7164cc28894Sstefano_zampini ierr = MatSetType(*C,deft);CHKERRQ(ierr); 7174cc28894Sstefano_zampini } 7184cc28894Sstefano_zampini (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_HYPRE; 7194cc28894Sstefano_zampini ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 7204cc28894Sstefano_zampini } 7214cc28894Sstefano_zampini ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 7224cc28894Sstefano_zampini ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 723613e5ff0Sstefano_zampini ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 724c1a070e6SStefano Zampini PetscFunctionReturn(0); 725c1a070e6SStefano Zampini } 726c1a070e6SStefano Zampini 7274cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_HYPRE_HYPRE(Mat A,Mat P,Mat C) 7284cc28894Sstefano_zampini { 7294cc28894Sstefano_zampini Mat B; 7304cc28894Sstefano_zampini hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr; 7314cc28894Sstefano_zampini Mat_HYPRE *hA,*hP; 7324cc28894Sstefano_zampini PetscBool ishypre; 7334cc28894Sstefano_zampini HYPRE_Int type; 7344cc28894Sstefano_zampini PetscErrorCode ierr; 7354cc28894Sstefano_zampini 7364cc28894Sstefano_zampini PetscFunctionBegin; 7374cc28894Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr); 7384cc28894Sstefano_zampini if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)P),PETSC_ERR_USER,"P should be of type %s",MATHYPRE); 7394cc28894Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr); 7404cc28894Sstefano_zampini if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE); 7414cc28894Sstefano_zampini hA = (Mat_HYPRE*)A->data; 7424cc28894Sstefano_zampini hP = (Mat_HYPRE*)P->data; 7434cc28894Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 7444cc28894Sstefano_zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 7454cc28894Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type)); 7464cc28894Sstefano_zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 7474cc28894Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr)); 7484cc28894Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr)); 7494cc28894Sstefano_zampini ierr = MatHYPRE_ParCSR_RAP(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr);CHKERRQ(ierr); 7504cc28894Sstefano_zampini ierr = MatCreateFromParCSR(ptapparcsr,MATHYPRE,PETSC_OWN_POINTER,&B);CHKERRQ(ierr); 7514cc28894Sstefano_zampini ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr); 7524cc28894Sstefano_zampini PetscFunctionReturn(0); 7534cc28894Sstefano_zampini } 7544cc28894Sstefano_zampini 755cd8bc7baSStefano Zampini static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 756cd8bc7baSStefano Zampini { 757cd8bc7baSStefano Zampini PetscErrorCode ierr; 758cd8bc7baSStefano Zampini 759cd8bc7baSStefano Zampini PetscFunctionBegin; 7604cc28894Sstefano_zampini if (scall == MAT_INITIAL_MATRIX) { 7614cc28894Sstefano_zampini ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 7624cc28894Sstefano_zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr); 7634cc28894Sstefano_zampini ierr = MatSetType(*C,MATHYPRE);CHKERRQ(ierr); 7644cc28894Sstefano_zampini (*C)->ops->ptapnumeric = MatPtAPNumeric_HYPRE_HYPRE; 7654cc28894Sstefano_zampini ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 7664cc28894Sstefano_zampini } 767613e5ff0Sstefano_zampini ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 7684cc28894Sstefano_zampini ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 769613e5ff0Sstefano_zampini ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 770cd8bc7baSStefano Zampini PetscFunctionReturn(0); 771cd8bc7baSStefano Zampini } 772cd8bc7baSStefano Zampini 773d501dc42Sstefano_zampini /* calls hypre_ParMatmul 774d501dc42Sstefano_zampini hypre_ParMatMul uses hypre_ParMatrixCreate with the communicator of hA 7753dad0653Sstefano_zampini hypre_ParMatrixCreate does not duplicate the communicator 7763dad0653Sstefano_zampini It looks like we don't need to have the diagonal entries 7773dad0653Sstefano_zampini ordered first in the rows of the diagonal part 7783dad0653Sstefano_zampini for boomerAMGBuildCoarseOperator to work */ 779d501dc42Sstefano_zampini static PetscErrorCode MatHYPRE_ParCSR_MatMatMult(hypre_ParCSRMatrix *hA, hypre_ParCSRMatrix *hB, hypre_ParCSRMatrix **hAB) 780d501dc42Sstefano_zampini { 781d501dc42Sstefano_zampini PetscFunctionBegin; 782d501dc42Sstefano_zampini PetscStackPush("hypre_ParMatmul"); 783d501dc42Sstefano_zampini *hAB = hypre_ParMatmul(hA,hB); 784d501dc42Sstefano_zampini PetscStackPop; 785d501dc42Sstefano_zampini PetscFunctionReturn(0); 786d501dc42Sstefano_zampini } 787d501dc42Sstefano_zampini 7885e5acdf2Sstefano_zampini static PetscErrorCode MatMatMultNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C) 7895e5acdf2Sstefano_zampini { 7905e5acdf2Sstefano_zampini Mat D; 791d501dc42Sstefano_zampini hypre_ParCSRMatrix *hA,*hB,*hAB = NULL; 7925e5acdf2Sstefano_zampini PetscErrorCode ierr; 7935e5acdf2Sstefano_zampini 7945e5acdf2Sstefano_zampini PetscFunctionBegin; 7955e5acdf2Sstefano_zampini ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr); 7965e5acdf2Sstefano_zampini ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr); 797d501dc42Sstefano_zampini ierr = MatHYPRE_ParCSR_MatMatMult(hA,hB,&hAB);CHKERRQ(ierr); 7985e5acdf2Sstefano_zampini ierr = MatCreateFromParCSR(hAB,MATAIJ,PETSC_OWN_POINTER,&D);CHKERRQ(ierr); 7995e5acdf2Sstefano_zampini ierr = MatHeaderMerge(C,&D);CHKERRQ(ierr); 8005e5acdf2Sstefano_zampini ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr); 8015e5acdf2Sstefano_zampini ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr); 8025e5acdf2Sstefano_zampini PetscFunctionReturn(0); 8035e5acdf2Sstefano_zampini } 8045e5acdf2Sstefano_zampini 8055e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat B,PetscReal fill,Mat *C) 8065e5acdf2Sstefano_zampini { 8075e5acdf2Sstefano_zampini PetscErrorCode ierr; 80820e1dc0dSstefano_zampini 8095e5acdf2Sstefano_zampini PetscFunctionBegin; 8105e5acdf2Sstefano_zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr); 8115e5acdf2Sstefano_zampini ierr = MatSetType(*C,MATAIJ);CHKERRQ(ierr); 8125e5acdf2Sstefano_zampini (*C)->ops->matmultnumeric = MatMatMultNumeric_AIJ_AIJ_wHYPRE; 8135e5acdf2Sstefano_zampini PetscFunctionReturn(0); 8145e5acdf2Sstefano_zampini } 8155e5acdf2Sstefano_zampini 816d501dc42Sstefano_zampini static PetscErrorCode MatMatMultNumeric_HYPRE_HYPRE(Mat A,Mat B,Mat C) 817d501dc42Sstefano_zampini { 818d501dc42Sstefano_zampini Mat D; 819d501dc42Sstefano_zampini hypre_ParCSRMatrix *Aparcsr,*Bparcsr,*ABparcsr = NULL; 820d501dc42Sstefano_zampini Mat_HYPRE *hA,*hB; 821d501dc42Sstefano_zampini PetscBool ishypre; 822d501dc42Sstefano_zampini HYPRE_Int type; 823d501dc42Sstefano_zampini PetscErrorCode ierr; 824d501dc42Sstefano_zampini 825d501dc42Sstefano_zampini PetscFunctionBegin; 826d501dc42Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)B,MATHYPRE,&ishypre);CHKERRQ(ierr); 827d501dc42Sstefano_zampini if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_USER,"B should be of type %s",MATHYPRE); 828d501dc42Sstefano_zampini ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr); 829d501dc42Sstefano_zampini if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE); 830d501dc42Sstefano_zampini hA = (Mat_HYPRE*)A->data; 831d501dc42Sstefano_zampini hB = (Mat_HYPRE*)B->data; 832d501dc42Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 833d501dc42Sstefano_zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 834d501dc42Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hB->ij,&type)); 835d501dc42Sstefano_zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported"); 836d501dc42Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr)); 837d501dc42Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hB->ij,(void**)&Bparcsr)); 838d501dc42Sstefano_zampini ierr = MatHYPRE_ParCSR_MatMatMult(Aparcsr,Bparcsr,&ABparcsr);CHKERRQ(ierr); 839d501dc42Sstefano_zampini ierr = MatCreateFromParCSR(ABparcsr,MATHYPRE,PETSC_OWN_POINTER,&D);CHKERRQ(ierr); 840d501dc42Sstefano_zampini /* need to use HeaderReplace because HeaderMerge messes up with the communicator */ 841d501dc42Sstefano_zampini ierr = MatHeaderReplace(C,&D);CHKERRQ(ierr); 842d501dc42Sstefano_zampini C->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE; 843d501dc42Sstefano_zampini PetscFunctionReturn(0); 844d501dc42Sstefano_zampini } 845d501dc42Sstefano_zampini 846d501dc42Sstefano_zampini static PetscErrorCode MatMatMult_HYPRE_HYPRE(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 847d501dc42Sstefano_zampini { 848d501dc42Sstefano_zampini PetscErrorCode ierr; 849d501dc42Sstefano_zampini 850d501dc42Sstefano_zampini PetscFunctionBegin; 851d501dc42Sstefano_zampini if (scall == MAT_INITIAL_MATRIX) { 852d501dc42Sstefano_zampini ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 853d501dc42Sstefano_zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr); 854d501dc42Sstefano_zampini ierr = MatSetType(*C,MATHYPRE);CHKERRQ(ierr); 855d501dc42Sstefano_zampini (*C)->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE; 856d501dc42Sstefano_zampini ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 857d501dc42Sstefano_zampini } 858d501dc42Sstefano_zampini ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 859d501dc42Sstefano_zampini ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 860d501dc42Sstefano_zampini ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 861d501dc42Sstefano_zampini PetscFunctionReturn(0); 862d501dc42Sstefano_zampini } 863d501dc42Sstefano_zampini 8643dad0653Sstefano_zampini PETSC_INTERN PetscErrorCode MatTransposeMatMatMultNumeric_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,Mat D) 86520e1dc0dSstefano_zampini { 86620e1dc0dSstefano_zampini Mat E; 86720e1dc0dSstefano_zampini hypre_ParCSRMatrix *hA,*hB,*hC,*hABC; 86820e1dc0dSstefano_zampini PetscErrorCode ierr; 86920e1dc0dSstefano_zampini 87020e1dc0dSstefano_zampini PetscFunctionBegin; 87120e1dc0dSstefano_zampini ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr); 87220e1dc0dSstefano_zampini ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr); 87320e1dc0dSstefano_zampini ierr = MatAIJGetParCSR_Private(C,&hC);CHKERRQ(ierr); 87420e1dc0dSstefano_zampini ierr = MatHYPRE_ParCSR_RAP(hA,hB,hC,&hABC);CHKERRQ(ierr); 87520e1dc0dSstefano_zampini ierr = MatCreateFromParCSR(hABC,MATAIJ,PETSC_OWN_POINTER,&E);CHKERRQ(ierr); 87620e1dc0dSstefano_zampini ierr = MatHeaderMerge(D,&E);CHKERRQ(ierr); 87720e1dc0dSstefano_zampini ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr); 87820e1dc0dSstefano_zampini ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr); 87920e1dc0dSstefano_zampini ierr = MatAIJRestoreParCSR_Private(C,&hC);CHKERRQ(ierr); 88020e1dc0dSstefano_zampini PetscFunctionReturn(0); 88120e1dc0dSstefano_zampini } 88220e1dc0dSstefano_zampini 8833dad0653Sstefano_zampini PETSC_INTERN PetscErrorCode MatTransposeMatMatMultSymbolic_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,PetscReal fill,Mat *D) 88420e1dc0dSstefano_zampini { 88520e1dc0dSstefano_zampini PetscErrorCode ierr; 88620e1dc0dSstefano_zampini 88720e1dc0dSstefano_zampini PetscFunctionBegin; 88820e1dc0dSstefano_zampini ierr = MatCreate(PetscObjectComm((PetscObject)A),D);CHKERRQ(ierr); 88920e1dc0dSstefano_zampini ierr = MatSetType(*D,MATAIJ);CHKERRQ(ierr); 89020e1dc0dSstefano_zampini PetscFunctionReturn(0); 89120e1dc0dSstefano_zampini } 89220e1dc0dSstefano_zampini 893ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y) 89463c07aadSStefano Zampini { 89563c07aadSStefano Zampini PetscErrorCode ierr; 89663c07aadSStefano Zampini 89763c07aadSStefano Zampini PetscFunctionBegin; 898414bd5c3SStefano Zampini ierr = MatHYPRE_MultKernel_Private(A,1.0,x,0.0,y,PETSC_TRUE);CHKERRQ(ierr); 89963c07aadSStefano Zampini PetscFunctionReturn(0); 90063c07aadSStefano Zampini } 90163c07aadSStefano Zampini 902ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y) 90363c07aadSStefano Zampini { 90463c07aadSStefano Zampini PetscErrorCode ierr; 90563c07aadSStefano Zampini 90663c07aadSStefano Zampini PetscFunctionBegin; 907414bd5c3SStefano Zampini ierr = MatHYPRE_MultKernel_Private(A,1.0,x,0.0,y,PETSC_FALSE);CHKERRQ(ierr); 90863c07aadSStefano Zampini PetscFunctionReturn(0); 90963c07aadSStefano Zampini } 91063c07aadSStefano Zampini 911414bd5c3SStefano Zampini static PetscErrorCode MatMultAdd_HYPRE(Mat A, Vec x, Vec y, Vec z) 912414bd5c3SStefano Zampini { 913414bd5c3SStefano Zampini PetscErrorCode ierr; 914414bd5c3SStefano Zampini 915414bd5c3SStefano Zampini PetscFunctionBegin; 916414bd5c3SStefano Zampini if (y != z) { 917414bd5c3SStefano Zampini ierr = VecCopy(y,z);CHKERRQ(ierr); 918414bd5c3SStefano Zampini } 919414bd5c3SStefano Zampini ierr = MatHYPRE_MultKernel_Private(A,1.0,x,1.0,z,PETSC_FALSE);CHKERRQ(ierr); 920414bd5c3SStefano Zampini PetscFunctionReturn(0); 921414bd5c3SStefano Zampini } 922414bd5c3SStefano Zampini 923414bd5c3SStefano Zampini static PetscErrorCode MatMultTransposeAdd_HYPRE(Mat A, Vec x, Vec y, Vec z) 924414bd5c3SStefano Zampini { 925414bd5c3SStefano Zampini PetscErrorCode ierr; 926414bd5c3SStefano Zampini 927414bd5c3SStefano Zampini PetscFunctionBegin; 928414bd5c3SStefano Zampini if (y != z) { 929414bd5c3SStefano Zampini ierr = VecCopy(y,z);CHKERRQ(ierr); 930414bd5c3SStefano Zampini } 931414bd5c3SStefano Zampini ierr = MatHYPRE_MultKernel_Private(A,1.0,x,1.0,z,PETSC_TRUE);CHKERRQ(ierr); 932414bd5c3SStefano Zampini PetscFunctionReturn(0); 933414bd5c3SStefano Zampini } 934414bd5c3SStefano Zampini 935414bd5c3SStefano Zampini /* y = a * A * x + b * y or y = a * A^t * x + b * y depending on trans */ 936414bd5c3SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, PetscScalar a, Vec x, PetscScalar b, Vec y, PetscBool trans) 93763c07aadSStefano Zampini { 93863c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 93963c07aadSStefano Zampini hypre_ParCSRMatrix *parcsr; 94063c07aadSStefano Zampini hypre_ParVector *hx,*hy; 94163c07aadSStefano Zampini PetscScalar *ax,*ay,*sax,*say; 94263c07aadSStefano Zampini PetscErrorCode ierr; 94363c07aadSStefano Zampini 94463c07aadSStefano Zampini PetscFunctionBegin; 94563c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr)); 94663c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx)); 94763c07aadSStefano Zampini PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy)); 94863c07aadSStefano Zampini ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); 94963c07aadSStefano Zampini ierr = VecGetArray(y,&ay);CHKERRQ(ierr); 95063c07aadSStefano Zampini if (trans) { 95158968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,ay,say); 95258968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax); 953414bd5c3SStefano Zampini hypre_ParCSRMatrixMatvecT(a,parcsr,hy,b,hx); 95458968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,say,ay); 95558968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax); 95663c07aadSStefano Zampini } else { 95758968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax); 95858968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,ay,say); 959414bd5c3SStefano Zampini hypre_ParCSRMatrixMatvec(a,parcsr,hx,b,hy); 96058968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax); 96158968eb6SStefano Zampini VecHYPRE_ParVectorReplacePointer(hA->b,say,ay); 96263c07aadSStefano Zampini } 96363c07aadSStefano Zampini ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr); 96463c07aadSStefano Zampini ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr); 96563c07aadSStefano Zampini PetscFunctionReturn(0); 96663c07aadSStefano Zampini } 96763c07aadSStefano Zampini 968ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A) 96963c07aadSStefano Zampini { 97063c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 97163c07aadSStefano Zampini PetscErrorCode ierr; 97263c07aadSStefano Zampini 97363c07aadSStefano Zampini PetscFunctionBegin; 97463c07aadSStefano Zampini if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x)); 97563c07aadSStefano Zampini if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b)); 976978814f1SStefano Zampini if (hA->ij) { 977978814f1SStefano Zampini if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL; 978978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij)); 979978814f1SStefano Zampini } 98063c07aadSStefano Zampini if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); } 981c69f721fSFande Kong 982c69f721fSFande Kong ierr = MatStashDestroy_Private(&A->stash);CHKERRQ(ierr); 983c69f721fSFande Kong 984c69f721fSFande Kong ierr = PetscFree(hA->array);CHKERRQ(ierr); 985c69f721fSFande Kong 98663c07aadSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr); 9872df22349SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr); 988c1a070e6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr); 989c1a070e6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr); 990d975228cSstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPRESetPreallocation_C",NULL);CHKERRQ(ierr); 991dd9c0a25Sstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPREGetParCSR_C",NULL);CHKERRQ(ierr); 99275d48cdbSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_is_hypre_C",NULL);CHKERRQ(ierr); 99363c07aadSStefano Zampini ierr = PetscFree(A->data);CHKERRQ(ierr); 99463c07aadSStefano Zampini PetscFunctionReturn(0); 99563c07aadSStefano Zampini } 99663c07aadSStefano Zampini 997ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A) 99863c07aadSStefano Zampini { 9994ec6421dSstefano_zampini PetscErrorCode ierr; 10004ec6421dSstefano_zampini 10014ec6421dSstefano_zampini PetscFunctionBegin; 10024ec6421dSstefano_zampini ierr = MatHYPRESetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr); 10034ec6421dSstefano_zampini PetscFunctionReturn(0); 10044ec6421dSstefano_zampini } 10054ec6421dSstefano_zampini 10064ec6421dSstefano_zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode) 10074ec6421dSstefano_zampini { 100863c07aadSStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 100963c07aadSStefano Zampini Vec x,b; 1010c69f721fSFande Kong PetscMPIInt n; 1011c69f721fSFande Kong PetscInt i,j,rstart,ncols,flg; 1012c69f721fSFande Kong PetscInt *row,*col; 1013c69f721fSFande Kong PetscScalar *val; 101463c07aadSStefano Zampini PetscErrorCode ierr; 101563c07aadSStefano Zampini 101663c07aadSStefano Zampini PetscFunctionBegin; 10174ec6421dSstefano_zampini if (mode == MAT_FLUSH_ASSEMBLY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_FLUSH_ASSEMBLY currently not supported with MATHYPRE"); 1018c69f721fSFande Kong 1019c69f721fSFande Kong if (!A->nooffprocentries) { 1020c69f721fSFande Kong while (1) { 1021c69f721fSFande Kong ierr = MatStashScatterGetMesg_Private(&A->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr); 1022c69f721fSFande Kong if (!flg) break; 1023c69f721fSFande Kong 1024c69f721fSFande Kong for (i=0; i<n; ) { 1025c69f721fSFande Kong /* Now identify the consecutive vals belonging to the same row */ 1026c69f721fSFande Kong for (j=i,rstart=row[j]; j<n; j++) { 1027c69f721fSFande Kong if (row[j] != rstart) break; 1028c69f721fSFande Kong } 1029c69f721fSFande Kong if (j < n) ncols = j-i; 1030c69f721fSFande Kong else ncols = n-i; 1031c69f721fSFande Kong /* Now assemble all these values with a single function call */ 1032c69f721fSFande Kong ierr = MatSetValues_HYPRE(A,1,row+i,ncols,col+i,val+i,A->insertmode);CHKERRQ(ierr); 1033c69f721fSFande Kong 1034c69f721fSFande Kong i = j; 1035c69f721fSFande Kong } 1036c69f721fSFande Kong } 1037c69f721fSFande Kong ierr = MatStashScatterEnd_Private(&A->stash);CHKERRQ(ierr); 1038c69f721fSFande Kong } 1039c69f721fSFande Kong 10404ec6421dSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij)); 1041af1cf968SStefano Zampini /* The assembly routine destroys the aux_matrix, we recreate it here by calling HYPRE_IJMatrixInitalize */ 1042af1cf968SStefano Zampini { 1043af1cf968SStefano Zampini hypre_AuxParCSRMatrix *aux_matrix; 1044af1cf968SStefano Zampini 1045af1cf968SStefano Zampini /* call destroy just to make sure we do not leak anything */ 1046af1cf968SStefano Zampini aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij); 1047af1cf968SStefano Zampini PetscStackCallStandard(hypre_AuxParCSRMatrixDestroy,(aux_matrix)); 1048af1cf968SStefano Zampini hypre_IJMatrixTranslator(hA->ij) = NULL; 1049af1cf968SStefano Zampini 1050af1cf968SStefano Zampini /* Initialize with assembled flag -> it only recreates the aux_par_matrix */ 1051af1cf968SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij)); 1052af1cf968SStefano Zampini aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij); 1053af1cf968SStefano Zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1; /* see comment in MatHYPRESetPreallocation_HYPRE */ 1054af1cf968SStefano Zampini PetscStackCallStandard(hypre_AuxParCSRMatrixInitialize,(aux_matrix)); 1055af1cf968SStefano Zampini } 105663c07aadSStefano Zampini if (hA->x) PetscFunctionReturn(0); 105763c07aadSStefano Zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 105863c07aadSStefano Zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 105963c07aadSStefano Zampini ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr); 106063c07aadSStefano Zampini ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr); 106163c07aadSStefano Zampini ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr); 106263c07aadSStefano Zampini ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr); 106363c07aadSStefano Zampini ierr = VecDestroy(&x);CHKERRQ(ierr); 106463c07aadSStefano Zampini ierr = VecDestroy(&b);CHKERRQ(ierr); 106563c07aadSStefano Zampini PetscFunctionReturn(0); 106663c07aadSStefano Zampini } 106763c07aadSStefano Zampini 1068c69f721fSFande Kong static PetscErrorCode MatGetArray_HYPRE(Mat A, PetscInt size, void **array) 1069c69f721fSFande Kong { 1070c69f721fSFande Kong Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 1071c69f721fSFande Kong PetscErrorCode ierr; 1072c69f721fSFande Kong 1073c69f721fSFande Kong PetscFunctionBegin; 1074c69f721fSFande Kong if (!hA->available) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Temporary space is in use"); 1075c69f721fSFande Kong 1076c69f721fSFande Kong if (hA->size >= size) *array = hA->array; 1077c69f721fSFande Kong else { 1078c69f721fSFande Kong ierr = PetscFree(hA->array);CHKERRQ(ierr); 1079c69f721fSFande Kong hA->size = size; 1080c69f721fSFande Kong ierr = PetscMalloc(hA->size,&hA->array);CHKERRQ(ierr); 1081c69f721fSFande Kong *array = hA->array; 1082c69f721fSFande Kong } 1083c69f721fSFande Kong 1084c69f721fSFande Kong hA->available = PETSC_FALSE; 1085c69f721fSFande Kong PetscFunctionReturn(0); 1086c69f721fSFande Kong } 1087c69f721fSFande Kong 1088708542d2SFande Kong static PetscErrorCode MatRestoreArray_HYPRE(Mat A, void **array) 1089c69f721fSFande Kong { 1090c69f721fSFande Kong Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 1091c69f721fSFande Kong 1092c69f721fSFande Kong PetscFunctionBegin; 1093c69f721fSFande Kong *array = NULL; 1094c69f721fSFande Kong hA->available = PETSC_TRUE; 1095c69f721fSFande Kong PetscFunctionReturn(0); 1096c69f721fSFande Kong } 1097c69f721fSFande Kong 1098d975228cSstefano_zampini PetscErrorCode MatSetValues_HYPRE(Mat A, PetscInt nr, const PetscInt rows[], PetscInt nc, const PetscInt cols[], const PetscScalar v[], InsertMode ins) 1099d975228cSstefano_zampini { 1100d975228cSstefano_zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 1101d975228cSstefano_zampini PetscScalar *vals = (PetscScalar *)v; 1102c69f721fSFande Kong PetscScalar *sscr; 1103c69f721fSFande Kong PetscInt *cscr[2]; 1104c69f721fSFande Kong PetscInt i,nzc; 110508defe43SFande Kong void *array = NULL; 1106d975228cSstefano_zampini PetscErrorCode ierr; 1107d975228cSstefano_zampini 1108d975228cSstefano_zampini PetscFunctionBegin; 1109c69f721fSFande Kong ierr = MatGetArray_HYPRE(A,sizeof(PetscInt)*(2*nc)+sizeof(PetscScalar)*nc*nr,&array);CHKERRQ(ierr); 1110c69f721fSFande Kong cscr[0] = (PetscInt*)array; 1111c69f721fSFande Kong cscr[1] = ((PetscInt*)array)+nc; 1112c69f721fSFande Kong sscr = (PetscScalar*)(((PetscInt*)array)+nc*2); 1113d975228cSstefano_zampini for (i=0,nzc=0;i<nc;i++) { 1114d975228cSstefano_zampini if (cols[i] >= 0) { 1115d975228cSstefano_zampini cscr[0][nzc ] = cols[i]; 1116d975228cSstefano_zampini cscr[1][nzc++] = i; 1117d975228cSstefano_zampini } 1118d975228cSstefano_zampini } 1119c69f721fSFande Kong if (!nzc) { 1120708542d2SFande Kong ierr = MatRestoreArray_HYPRE(A,&array);CHKERRQ(ierr); 1121c69f721fSFande Kong PetscFunctionReturn(0); 1122c69f721fSFande Kong } 1123d975228cSstefano_zampini 1124d975228cSstefano_zampini if (ins == ADD_VALUES) { 1125d975228cSstefano_zampini for (i=0;i<nr;i++) { 1126e3977e59Sstefano_zampini if (rows[i] >= 0 && nzc) { 1127d975228cSstefano_zampini PetscInt j; 1128d975228cSstefano_zampini for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]]; 112908defe43SFande Kong PetscStackCallStandard(HYPRE_IJMatrixAddToValues,(hA->ij,1,(HYPRE_Int*)&nzc,(HYPRE_Int*)(rows+i),(HYPRE_Int*)cscr[0],sscr)); 1130d975228cSstefano_zampini } 1131d975228cSstefano_zampini vals += nc; 1132d975228cSstefano_zampini } 1133d975228cSstefano_zampini } else { /* INSERT_VALUES */ 1134c69f721fSFande Kong 1135d975228cSstefano_zampini PetscInt rst,ren; 1136d975228cSstefano_zampini ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr); 1137c69f721fSFande Kong 1138d975228cSstefano_zampini for (i=0;i<nr;i++) { 1139e3977e59Sstefano_zampini if (rows[i] >= 0 && nzc) { 1140d975228cSstefano_zampini PetscInt j; 1141d975228cSstefano_zampini for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]]; 1142c69f721fSFande Kong /* nonlocal values */ 1143c69f721fSFande Kong if (rows[i] < rst || rows[i] >= ren) { ierr = MatStashValuesRow_Private(&A->stash,rows[i],nzc,cscr[0],sscr,PETSC_FALSE);CHKERRQ(ierr); } 1144c69f721fSFande Kong /* local values */ 114508defe43SFande Kong else PetscStackCallStandard(HYPRE_IJMatrixSetValues,(hA->ij,1,(HYPRE_Int*)&nzc,(HYPRE_Int*)(rows+i),(HYPRE_Int*)cscr[0],sscr)); 1146d975228cSstefano_zampini } 1147d975228cSstefano_zampini vals += nc; 1148d975228cSstefano_zampini } 1149d975228cSstefano_zampini } 1150c69f721fSFande Kong 1151708542d2SFande Kong ierr = MatRestoreArray_HYPRE(A,&array);CHKERRQ(ierr); 1152d975228cSstefano_zampini PetscFunctionReturn(0); 1153d975228cSstefano_zampini } 1154d975228cSstefano_zampini 1155d975228cSstefano_zampini static PetscErrorCode MatHYPRESetPreallocation_HYPRE(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[]) 1156d975228cSstefano_zampini { 1157d975228cSstefano_zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 11587d968826Sstefano_zampini HYPRE_Int *hdnnz,*honnz; 115906a29025Sstefano_zampini PetscInt i,rs,re,cs,ce,bs; 1160d975228cSstefano_zampini PetscMPIInt size; 1161d975228cSstefano_zampini PetscErrorCode ierr; 1162d975228cSstefano_zampini 1163d975228cSstefano_zampini PetscFunctionBegin; 116406a29025Sstefano_zampini ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1165d975228cSstefano_zampini ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr); 1166d975228cSstefano_zampini ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr); 1167d975228cSstefano_zampini rs = A->rmap->rstart; 1168d975228cSstefano_zampini re = A->rmap->rend; 1169d975228cSstefano_zampini cs = A->cmap->rstart; 1170d975228cSstefano_zampini ce = A->cmap->rend; 1171d975228cSstefano_zampini if (!hA->ij) { 1172d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rs,re-1,cs,ce-1,&hA->ij)); 1173d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); 1174d975228cSstefano_zampini } else { 11757d968826Sstefano_zampini HYPRE_Int hrs,hre,hcs,hce; 1176d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetLocalRange,(hA->ij,&hrs,&hre,&hcs,&hce)); 1177d975228cSstefano_zampini if (hre-hrs+1 != re -rs) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent local rows: IJMatrix [%D,%D), PETSc [%D,%d)",hrs,hre+1,rs,re); 1178d975228cSstefano_zampini if (hce-hcs+1 != ce -cs) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent local cols: IJMatrix [%D,%D), PETSc [%D,%d)",hcs,hce+1,cs,ce); 1179d975228cSstefano_zampini } 118006a29025Sstefano_zampini if (dnz == PETSC_DEFAULT || dnz == PETSC_DECIDE) dnz = 10*bs; 118106a29025Sstefano_zampini if (onz == PETSC_DEFAULT || onz == PETSC_DECIDE) onz = 10*bs; 118206a29025Sstefano_zampini 1183d975228cSstefano_zampini if (!dnnz) { 1184d975228cSstefano_zampini ierr = PetscMalloc1(A->rmap->n,&hdnnz);CHKERRQ(ierr); 1185d975228cSstefano_zampini for (i=0;i<A->rmap->n;i++) hdnnz[i] = dnz; 1186d975228cSstefano_zampini } else { 11877d968826Sstefano_zampini hdnnz = (HYPRE_Int*)dnnz; 1188d975228cSstefano_zampini } 1189d975228cSstefano_zampini ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 1190d975228cSstefano_zampini if (size > 1) { 1191ddbeb582SStefano Zampini hypre_AuxParCSRMatrix *aux_matrix; 1192d975228cSstefano_zampini if (!onnz) { 1193d975228cSstefano_zampini ierr = PetscMalloc1(A->rmap->n,&honnz);CHKERRQ(ierr); 1194d975228cSstefano_zampini for (i=0;i<A->rmap->n;i++) honnz[i] = onz; 1195d975228cSstefano_zampini } else { 11967d968826Sstefano_zampini honnz = (HYPRE_Int*)onnz; 1197d975228cSstefano_zampini } 1198ddbeb582SStefano Zampini /* SetDiagOffdSizes sets hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0, since it seems 1199ddbeb582SStefano Zampini they assume the user will input the entire row values, properly sorted 1200ddbeb582SStefano Zampini In PETSc, we don't make such an assumption, and we instead set this flag to 1 1201ddbeb582SStefano Zampini Also, to avoid possible memory leaks, we destroy and recreate the translator 1202ddbeb582SStefano Zampini This has to be done here, as HYPRE_IJMatrixInitialize will properly initialize 1203ddbeb582SStefano Zampini the IJ matrix for us */ 1204ddbeb582SStefano Zampini aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij); 1205ddbeb582SStefano Zampini hypre_AuxParCSRMatrixDestroy(aux_matrix); 1206ddbeb582SStefano Zampini hypre_IJMatrixTranslator(hA->ij) = NULL; 1207d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(hA->ij,hdnnz,honnz)); 1208ddbeb582SStefano Zampini aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij); 1209ddbeb582SStefano Zampini hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1; 1210d975228cSstefano_zampini } else { 1211d975228cSstefano_zampini honnz = NULL; 1212d975228cSstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixSetRowSizes,(hA->ij,hdnnz)); 1213d975228cSstefano_zampini } 1214ddbeb582SStefano Zampini 1215af1cf968SStefano Zampini /* reset assembled flag and call the initialize method */ 1216af1cf968SStefano Zampini hypre_IJMatrixAssembleFlag(hA->ij) = 0; 1217ddbeb582SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij)); 1218d975228cSstefano_zampini if (!dnnz) { 1219d975228cSstefano_zampini ierr = PetscFree(hdnnz);CHKERRQ(ierr); 1220d975228cSstefano_zampini } 1221d975228cSstefano_zampini if (!onnz && honnz) { 1222d975228cSstefano_zampini ierr = PetscFree(honnz);CHKERRQ(ierr); 1223d975228cSstefano_zampini } 1224af1cf968SStefano Zampini 1225af1cf968SStefano Zampini /* Match AIJ logic */ 122606a29025Sstefano_zampini A->preallocated = PETSC_TRUE; 1227af1cf968SStefano Zampini A->assembled = PETSC_FALSE; 1228d975228cSstefano_zampini PetscFunctionReturn(0); 1229d975228cSstefano_zampini } 1230d975228cSstefano_zampini 1231d975228cSstefano_zampini /*@C 1232d975228cSstefano_zampini MatHYPRESetPreallocation - Preallocates memory for a sparse parallel matrix in HYPRE IJ format 1233d975228cSstefano_zampini 1234d975228cSstefano_zampini Collective on Mat 1235d975228cSstefano_zampini 1236d975228cSstefano_zampini Input Parameters: 1237d975228cSstefano_zampini + A - the matrix 1238d975228cSstefano_zampini . dnz - number of nonzeros per row in DIAGONAL portion of local submatrix 1239d975228cSstefano_zampini (same value is used for all local rows) 1240d975228cSstefano_zampini . dnnz - array containing the number of nonzeros in the various rows of the 1241d975228cSstefano_zampini DIAGONAL portion of the local submatrix (possibly different for each row) 1242d975228cSstefano_zampini or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure. 1243d975228cSstefano_zampini The size of this array is equal to the number of local rows, i.e 'm'. 1244d975228cSstefano_zampini For matrices that will be factored, you must leave room for (and set) 1245d975228cSstefano_zampini the diagonal entry even if it is zero. 1246d975228cSstefano_zampini . onz - number of nonzeros per row in the OFF-DIAGONAL portion of local 1247d975228cSstefano_zampini submatrix (same value is used for all local rows). 1248d975228cSstefano_zampini - onnz - array containing the number of nonzeros in the various rows of the 1249d975228cSstefano_zampini OFF-DIAGONAL portion of the local submatrix (possibly different for 1250d975228cSstefano_zampini each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero 1251d975228cSstefano_zampini structure. The size of this array is equal to the number 1252d975228cSstefano_zampini of local rows, i.e 'm'. 1253d975228cSstefano_zampini 125495452b02SPatrick Sanan Notes: 125595452b02SPatrick Sanan If the *nnz parameter is given then the *nz parameter is ignored; for sequential matrices, onz and onnz are ignored. 1256d975228cSstefano_zampini 1257d975228cSstefano_zampini Level: intermediate 1258d975228cSstefano_zampini 1259d975228cSstefano_zampini .keywords: matrix, aij, compressed row, sparse, parallel 1260d975228cSstefano_zampini 1261af1cf968SStefano Zampini .seealso: MatCreate(), MatMPIAIJSetPreallocation(), MATHYPRE 1262d975228cSstefano_zampini @*/ 1263d975228cSstefano_zampini PetscErrorCode MatHYPRESetPreallocation(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[]) 1264d975228cSstefano_zampini { 1265d975228cSstefano_zampini PetscErrorCode ierr; 1266d975228cSstefano_zampini 1267d975228cSstefano_zampini PetscFunctionBegin; 1268d975228cSstefano_zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1269d975228cSstefano_zampini PetscValidType(A,1); 1270d975228cSstefano_zampini ierr = PetscTryMethod(A,"MatHYPRESetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(A,dnz,dnnz,onz,onnz));CHKERRQ(ierr); 1271d975228cSstefano_zampini PetscFunctionReturn(0); 1272d975228cSstefano_zampini } 1273d975228cSstefano_zampini 1274225daaf8SStefano Zampini /* 1275225daaf8SStefano Zampini MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix 1276225daaf8SStefano Zampini 1277225daaf8SStefano Zampini Collective 1278225daaf8SStefano Zampini 1279225daaf8SStefano Zampini Input Parameters: 128045b8d346SStefano Zampini + parcsr - the pointer to the hypre_ParCSRMatrix 1281bb4689ddSStefano Zampini . mtype - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported. 1282225daaf8SStefano Zampini - copymode - PETSc copying options 1283225daaf8SStefano Zampini 1284225daaf8SStefano Zampini Output Parameter: 1285225daaf8SStefano Zampini . A - the matrix 1286225daaf8SStefano Zampini 1287225daaf8SStefano Zampini Level: intermediate 1288225daaf8SStefano Zampini 1289225daaf8SStefano Zampini .seealso: MatHYPRE, PetscCopyMode 1290225daaf8SStefano Zampini */ 129145b8d346SStefano Zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(hypre_ParCSRMatrix *parcsr, MatType mtype, PetscCopyMode copymode, Mat* A) 1292978814f1SStefano Zampini { 1293225daaf8SStefano Zampini Mat T; 1294978814f1SStefano Zampini Mat_HYPRE *hA; 1295978814f1SStefano Zampini MPI_Comm comm; 1296978814f1SStefano Zampini PetscInt rstart,rend,cstart,cend,M,N; 1297bb4689ddSStefano Zampini PetscBool isseqaij,ismpiaij,isaij,ishyp,isis; 1298978814f1SStefano Zampini PetscErrorCode ierr; 1299978814f1SStefano Zampini 1300978814f1SStefano Zampini PetscFunctionBegin; 1301978814f1SStefano Zampini comm = hypre_ParCSRMatrixComm(parcsr); 1302225daaf8SStefano Zampini ierr = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 1303225daaf8SStefano Zampini ierr = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr); 1304225daaf8SStefano Zampini ierr = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr); 1305225daaf8SStefano Zampini ierr = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr); 13068cfe8d00SStefano Zampini ierr = PetscStrcmp(mtype,MATIS,&isis);CHKERRQ(ierr); 1307225daaf8SStefano Zampini isaij = (PetscBool)(isseqaij || ismpiaij || isaij); 1308bb4689ddSStefano Zampini if (!isaij && !ishyp && !isis) SETERRQ6(comm,PETSC_ERR_SUP,"Unsupported MatType %s! Supported types are %s, %s, %s, %s, and %s",mtype,MATAIJ,MATSEQAIJ,MATMPIAIJ,MATIS,MATHYPRE); 1309978814f1SStefano Zampini /* access ParCSRMatrix */ 1310978814f1SStefano Zampini rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr); 1311978814f1SStefano Zampini rend = hypre_ParCSRMatrixLastRowIndex(parcsr); 1312978814f1SStefano Zampini cstart = hypre_ParCSRMatrixFirstColDiag(parcsr); 1313978814f1SStefano Zampini cend = hypre_ParCSRMatrixLastColDiag(parcsr); 1314978814f1SStefano Zampini M = hypre_ParCSRMatrixGlobalNumRows(parcsr); 1315978814f1SStefano Zampini N = hypre_ParCSRMatrixGlobalNumCols(parcsr); 1316978814f1SStefano Zampini 1317fa92c42cSstefano_zampini /* fix for empty local rows/columns */ 1318fa92c42cSstefano_zampini if (rend < rstart) rend = rstart; 1319fa92c42cSstefano_zampini if (cend < cstart) cend = cstart; 1320fa92c42cSstefano_zampini 1321e6471dc9SStefano Zampini /* PETSc convention */ 1322e6471dc9SStefano Zampini rend++; 1323e6471dc9SStefano Zampini cend++; 1324e6471dc9SStefano Zampini rend = PetscMin(rend,M); 1325e6471dc9SStefano Zampini cend = PetscMin(cend,N); 1326e6471dc9SStefano Zampini 1327978814f1SStefano Zampini /* create PETSc matrix with MatHYPRE */ 1328225daaf8SStefano Zampini ierr = MatCreate(comm,&T);CHKERRQ(ierr); 1329e6471dc9SStefano Zampini ierr = MatSetSizes(T,rend-rstart,cend-cstart,M,N);CHKERRQ(ierr); 1330225daaf8SStefano Zampini ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr); 1331225daaf8SStefano Zampini hA = (Mat_HYPRE*)(T->data); 1332978814f1SStefano Zampini 1333978814f1SStefano Zampini /* create HYPRE_IJMatrix */ 1334978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij)); 133545b8d346SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR)); 133645b8d346SStefano Zampini 133745b8d346SStefano Zampini /* create new ParCSR object if needed */ 133845b8d346SStefano Zampini if (ishyp && copymode == PETSC_COPY_VALUES) { 133945b8d346SStefano Zampini hypre_ParCSRMatrix *new_parcsr; 134045b8d346SStefano Zampini hypre_CSRMatrix *hdiag,*hoffd,*ndiag,*noffd; 134145b8d346SStefano Zampini 134245b8d346SStefano Zampini new_parcsr = hypre_ParCSRMatrixCompleteClone(parcsr); 134345b8d346SStefano Zampini hdiag = hypre_ParCSRMatrixDiag(parcsr); 134445b8d346SStefano Zampini hoffd = hypre_ParCSRMatrixOffd(parcsr); 134545b8d346SStefano Zampini ndiag = hypre_ParCSRMatrixDiag(new_parcsr); 134645b8d346SStefano Zampini noffd = hypre_ParCSRMatrixOffd(new_parcsr); 134745b8d346SStefano Zampini ierr = PetscMemcpy(hypre_CSRMatrixData(ndiag),hypre_CSRMatrixData(hdiag),hypre_CSRMatrixNumNonzeros(hdiag)*sizeof(PetscScalar));CHKERRQ(ierr); 134845b8d346SStefano Zampini ierr = PetscMemcpy(hypre_CSRMatrixData(noffd),hypre_CSRMatrixData(hoffd),hypre_CSRMatrixNumNonzeros(hoffd)*sizeof(PetscScalar));CHKERRQ(ierr); 134945b8d346SStefano Zampini parcsr = new_parcsr; 135045b8d346SStefano Zampini copymode = PETSC_OWN_POINTER; 135145b8d346SStefano Zampini } 1352978814f1SStefano Zampini 1353978814f1SStefano Zampini /* set ParCSR object */ 1354978814f1SStefano Zampini hypre_IJMatrixObject(hA->ij) = parcsr; 13554ec6421dSstefano_zampini T->preallocated = PETSC_TRUE; 1356978814f1SStefano Zampini 1357978814f1SStefano Zampini /* set assembled flag */ 1358978814f1SStefano Zampini hypre_IJMatrixAssembleFlag(hA->ij) = 1; 1359978814f1SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij)); 1360225daaf8SStefano Zampini if (ishyp) { 13616d2a658fSstefano_zampini PetscMPIInt myid = 0; 13626d2a658fSstefano_zampini 13636d2a658fSstefano_zampini /* make sure we always have row_starts and col_starts available */ 13646d2a658fSstefano_zampini if (HYPRE_AssumedPartitionCheck()) { 13656d2a658fSstefano_zampini ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr); 13666d2a658fSstefano_zampini } 13676d2a658fSstefano_zampini if (!hypre_ParCSRMatrixOwnsColStarts(parcsr)) { 13686d2a658fSstefano_zampini PetscLayout map; 13696d2a658fSstefano_zampini 13706d2a658fSstefano_zampini ierr = MatGetLayouts(T,NULL,&map);CHKERRQ(ierr); 13716d2a658fSstefano_zampini ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 13721c92d2d0Sstefano_zampini hypre_ParCSRMatrixColStarts(parcsr) = (HYPRE_Int*)(map->range + myid); 13736d2a658fSstefano_zampini } 13746d2a658fSstefano_zampini if (!hypre_ParCSRMatrixOwnsRowStarts(parcsr)) { 13756d2a658fSstefano_zampini PetscLayout map; 13766d2a658fSstefano_zampini 13776d2a658fSstefano_zampini ierr = MatGetLayouts(T,&map,NULL);CHKERRQ(ierr); 13786d2a658fSstefano_zampini ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 13791c92d2d0Sstefano_zampini hypre_ParCSRMatrixRowStarts(parcsr) = (HYPRE_Int*)(map->range + myid); 13806d2a658fSstefano_zampini } 1381978814f1SStefano Zampini /* prevent from freeing the pointer */ 1382978814f1SStefano Zampini if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE; 1383225daaf8SStefano Zampini *A = T; 1384978814f1SStefano Zampini ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1385978814f1SStefano Zampini ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1386bb4689ddSStefano Zampini } else if (isaij) { 1387bb4689ddSStefano Zampini if (copymode != PETSC_OWN_POINTER) { 1388225daaf8SStefano Zampini /* prevent from freeing the pointer */ 1389225daaf8SStefano Zampini hA->inner_free = PETSC_FALSE; 1390225daaf8SStefano Zampini ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr); 1391225daaf8SStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 1392225daaf8SStefano Zampini } else { /* AIJ return type with PETSC_OWN_POINTER */ 1393225daaf8SStefano Zampini ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr); 1394225daaf8SStefano Zampini *A = T; 1395225daaf8SStefano Zampini } 1396bb4689ddSStefano Zampini } else if (isis) { 13978cfe8d00SStefano Zampini ierr = MatConvert_HYPRE_IS(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr); 13988cfe8d00SStefano Zampini if (copymode != PETSC_OWN_POINTER) hA->inner_free = PETSC_FALSE; 1399bb4689ddSStefano Zampini ierr = MatDestroy(&T);CHKERRQ(ierr); 1400bb4689ddSStefano Zampini } 1401978814f1SStefano Zampini PetscFunctionReturn(0); 1402978814f1SStefano Zampini } 1403978814f1SStefano Zampini 140468ec7858SStefano Zampini static PetscErrorCode MatHYPREGetParCSR_HYPRE(Mat A, hypre_ParCSRMatrix **parcsr) 1405dd9c0a25Sstefano_zampini { 1406dd9c0a25Sstefano_zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 1407dd9c0a25Sstefano_zampini HYPRE_Int type; 1408dd9c0a25Sstefano_zampini 1409dd9c0a25Sstefano_zampini PetscFunctionBegin; 1410dd9c0a25Sstefano_zampini if (!hA->ij) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"HYPRE_IJMatrix not present"); 1411dd9c0a25Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type)); 1412dd9c0a25Sstefano_zampini if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"HYPRE_IJMatrix is not of type HYPRE_PARCSR"); 1413dd9c0a25Sstefano_zampini PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)parcsr)); 1414dd9c0a25Sstefano_zampini PetscFunctionReturn(0); 1415dd9c0a25Sstefano_zampini } 1416dd9c0a25Sstefano_zampini 1417dd9c0a25Sstefano_zampini /* 1418dd9c0a25Sstefano_zampini MatHYPREGetParCSR - Gets the pointer to the ParCSR matrix 1419dd9c0a25Sstefano_zampini 1420dd9c0a25Sstefano_zampini Not collective 1421dd9c0a25Sstefano_zampini 1422dd9c0a25Sstefano_zampini Input Parameters: 1423dd9c0a25Sstefano_zampini + A - the MATHYPRE object 1424dd9c0a25Sstefano_zampini 1425dd9c0a25Sstefano_zampini Output Parameter: 1426dd9c0a25Sstefano_zampini . parcsr - the pointer to the hypre_ParCSRMatrix 1427dd9c0a25Sstefano_zampini 1428dd9c0a25Sstefano_zampini Level: intermediate 1429dd9c0a25Sstefano_zampini 1430dd9c0a25Sstefano_zampini .seealso: MatHYPRE, PetscCopyMode 1431dd9c0a25Sstefano_zampini */ 1432dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR(Mat A, hypre_ParCSRMatrix **parcsr) 1433dd9c0a25Sstefano_zampini { 1434dd9c0a25Sstefano_zampini PetscErrorCode ierr; 1435dd9c0a25Sstefano_zampini 1436dd9c0a25Sstefano_zampini PetscFunctionBegin; 1437dd9c0a25Sstefano_zampini PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1438dd9c0a25Sstefano_zampini PetscValidType(A,1); 1439dd9c0a25Sstefano_zampini ierr = PetscUseMethod(A,"MatHYPREGetParCSR_C",(Mat,hypre_ParCSRMatrix**),(A,parcsr));CHKERRQ(ierr); 1440dd9c0a25Sstefano_zampini PetscFunctionReturn(0); 1441dd9c0a25Sstefano_zampini } 1442dd9c0a25Sstefano_zampini 144368ec7858SStefano Zampini static PetscErrorCode MatMissingDiagonal_HYPRE(Mat A, PetscBool *missing, PetscInt *dd) 144468ec7858SStefano Zampini { 144568ec7858SStefano Zampini hypre_ParCSRMatrix *parcsr; 144668ec7858SStefano Zampini hypre_CSRMatrix *ha; 144768ec7858SStefano Zampini PetscInt rst; 144868ec7858SStefano Zampini PetscErrorCode ierr; 144968ec7858SStefano Zampini 145068ec7858SStefano Zampini PetscFunctionBegin; 145168ec7858SStefano Zampini if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented with non-square diagonal blocks"); 145268ec7858SStefano Zampini ierr = MatGetOwnershipRange(A,&rst,NULL);CHKERRQ(ierr); 145368ec7858SStefano Zampini ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr); 145468ec7858SStefano Zampini if (missing) *missing = PETSC_FALSE; 145568ec7858SStefano Zampini if (dd) *dd = -1; 145668ec7858SStefano Zampini ha = hypre_ParCSRMatrixDiag(parcsr); 145768ec7858SStefano Zampini if (ha) { 145868299464SStefano Zampini PetscInt size,i; 145968299464SStefano Zampini HYPRE_Int *ii,*jj; 146068ec7858SStefano Zampini 146168ec7858SStefano Zampini size = hypre_CSRMatrixNumRows(ha); 146268ec7858SStefano Zampini ii = hypre_CSRMatrixI(ha); 146368ec7858SStefano Zampini jj = hypre_CSRMatrixJ(ha); 146468ec7858SStefano Zampini for (i = 0; i < size; i++) { 146568ec7858SStefano Zampini PetscInt j; 146668ec7858SStefano Zampini PetscBool found = PETSC_FALSE; 146768ec7858SStefano Zampini 146868ec7858SStefano Zampini for (j = ii[i]; j < ii[i+1] && !found; j++) 146968ec7858SStefano Zampini found = (jj[j] == i) ? PETSC_TRUE : PETSC_FALSE; 147068ec7858SStefano Zampini 147168ec7858SStefano Zampini if (!found) { 147268ec7858SStefano Zampini PetscInfo1(A,"Matrix is missing local diagonal entry %D\n",i); 147368ec7858SStefano Zampini if (missing) *missing = PETSC_TRUE; 147468ec7858SStefano Zampini if (dd) *dd = i+rst; 147568ec7858SStefano Zampini PetscFunctionReturn(0); 147668ec7858SStefano Zampini } 147768ec7858SStefano Zampini } 147868ec7858SStefano Zampini if (!size) { 147968ec7858SStefano Zampini PetscInfo(A,"Matrix has no diagonal entries therefore is missing diagonal\n"); 148068ec7858SStefano Zampini if (missing) *missing = PETSC_TRUE; 148168ec7858SStefano Zampini if (dd) *dd = rst; 148268ec7858SStefano Zampini } 148368ec7858SStefano Zampini } else { 148468ec7858SStefano Zampini PetscInfo(A,"Matrix has no diagonal entries therefore is missing diagonal\n"); 148568ec7858SStefano Zampini if (missing) *missing = PETSC_TRUE; 148668ec7858SStefano Zampini if (dd) *dd = rst; 148768ec7858SStefano Zampini } 148868ec7858SStefano Zampini PetscFunctionReturn(0); 148968ec7858SStefano Zampini } 149068ec7858SStefano Zampini 149168ec7858SStefano Zampini static PetscErrorCode MatScale_HYPRE(Mat A, PetscScalar s) 149268ec7858SStefano Zampini { 149368ec7858SStefano Zampini hypre_ParCSRMatrix *parcsr; 149468ec7858SStefano Zampini hypre_CSRMatrix *ha; 149568ec7858SStefano Zampini PetscErrorCode ierr; 149668ec7858SStefano Zampini 149768ec7858SStefano Zampini PetscFunctionBegin; 149868ec7858SStefano Zampini ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr); 149968ec7858SStefano Zampini /* diagonal part */ 150068ec7858SStefano Zampini ha = hypre_ParCSRMatrixDiag(parcsr); 150168ec7858SStefano Zampini if (ha) { 150268299464SStefano Zampini PetscInt size,i; 150368299464SStefano Zampini HYPRE_Int *ii; 150468ec7858SStefano Zampini PetscScalar *a; 150568ec7858SStefano Zampini 150668ec7858SStefano Zampini size = hypre_CSRMatrixNumRows(ha); 150768ec7858SStefano Zampini a = hypre_CSRMatrixData(ha); 150868ec7858SStefano Zampini ii = hypre_CSRMatrixI(ha); 150968ec7858SStefano Zampini for (i = 0; i < ii[size]; i++) a[i] *= s; 151068ec7858SStefano Zampini } 151168ec7858SStefano Zampini /* offdiagonal part */ 151268ec7858SStefano Zampini ha = hypre_ParCSRMatrixOffd(parcsr); 151368ec7858SStefano Zampini if (ha) { 151468299464SStefano Zampini PetscInt size,i; 151568299464SStefano Zampini HYPRE_Int *ii; 151668ec7858SStefano Zampini PetscScalar *a; 151768ec7858SStefano Zampini 151868ec7858SStefano Zampini size = hypre_CSRMatrixNumRows(ha); 151968ec7858SStefano Zampini a = hypre_CSRMatrixData(ha); 152068ec7858SStefano Zampini ii = hypre_CSRMatrixI(ha); 152168ec7858SStefano Zampini for (i = 0; i < ii[size]; i++) a[i] *= s; 152268ec7858SStefano Zampini } 152368ec7858SStefano Zampini PetscFunctionReturn(0); 152468ec7858SStefano Zampini } 152568ec7858SStefano Zampini 152668ec7858SStefano Zampini static PetscErrorCode MatZeroRowsColumns_HYPRE(Mat A, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b) 152768ec7858SStefano Zampini { 152868ec7858SStefano Zampini hypre_ParCSRMatrix *parcsr; 152968299464SStefano Zampini HYPRE_Int *lrows; 153068299464SStefano Zampini PetscInt rst,ren,i; 153168ec7858SStefano Zampini PetscErrorCode ierr; 153268ec7858SStefano Zampini 153368ec7858SStefano Zampini PetscFunctionBegin; 153468ec7858SStefano Zampini if (x || b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented"); 153568ec7858SStefano Zampini ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr); 153668ec7858SStefano Zampini ierr = PetscMalloc1(numRows,&lrows);CHKERRQ(ierr); 153768ec7858SStefano Zampini ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr); 153868ec7858SStefano Zampini for (i=0;i<numRows;i++) { 153968ec7858SStefano Zampini if (rows[i] < rst || rows[i] >= ren) 154068ec7858SStefano Zampini SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Non-local rows not yet supported"); 154168ec7858SStefano Zampini lrows[i] = rows[i] - rst; 154268ec7858SStefano Zampini } 154368ec7858SStefano Zampini PetscStackCallStandard(hypre_ParCSRMatrixEliminateRowsCols,(parcsr,numRows,lrows)); 154468ec7858SStefano Zampini ierr = PetscFree(lrows);CHKERRQ(ierr); 154568ec7858SStefano Zampini PetscFunctionReturn(0); 154668ec7858SStefano Zampini } 154768ec7858SStefano Zampini 1548c69f721fSFande Kong static PetscErrorCode MatZeroEntries_HYPRE_CSRMatrix(hypre_CSRMatrix *ha) 1549c69f721fSFande Kong { 1550c69f721fSFande Kong PetscErrorCode ierr; 1551c69f721fSFande Kong 1552c69f721fSFande Kong PetscFunctionBegin; 1553c69f721fSFande Kong if (ha) { 1554c69f721fSFande Kong HYPRE_Int *ii, size; 1555c69f721fSFande Kong HYPRE_Complex *a; 1556c69f721fSFande Kong 1557c69f721fSFande Kong size = hypre_CSRMatrixNumRows(ha); 1558c69f721fSFande Kong a = hypre_CSRMatrixData(ha); 1559c69f721fSFande Kong ii = hypre_CSRMatrixI(ha); 1560c69f721fSFande Kong 1561c69f721fSFande Kong if (a) { ierr = PetscMemzero(a,(ii[size])*sizeof(HYPRE_Complex));CHKERRQ(ierr); } 1562c69f721fSFande Kong } 1563c69f721fSFande Kong PetscFunctionReturn(0); 1564c69f721fSFande Kong } 1565c69f721fSFande Kong 1566c69f721fSFande Kong PetscErrorCode MatZeroEntries_HYPRE(Mat A) 1567c69f721fSFande Kong { 1568c69f721fSFande Kong hypre_ParCSRMatrix *parcsr; 1569c69f721fSFande Kong PetscErrorCode ierr; 1570c69f721fSFande Kong 1571c69f721fSFande Kong PetscFunctionBegin; 1572c69f721fSFande Kong ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr); 1573c69f721fSFande Kong /* diagonal part */ 1574c69f721fSFande Kong ierr = MatZeroEntries_HYPRE_CSRMatrix(hypre_ParCSRMatrixDiag(parcsr));CHKERRQ(ierr); 1575c69f721fSFande Kong /* off-diagonal part */ 1576c69f721fSFande Kong ierr = MatZeroEntries_HYPRE_CSRMatrix(hypre_ParCSRMatrixOffd(parcsr));CHKERRQ(ierr); 1577c69f721fSFande Kong PetscFunctionReturn(0); 1578c69f721fSFande Kong } 1579c69f721fSFande Kong 1580c69f721fSFande Kong static PetscErrorCode MatZeroRows_HYPRE_CSRMatrix(hypre_CSRMatrix *hA,PetscInt N,const PetscInt rows[],PetscScalar diag) 1581c69f721fSFande Kong { 1582c69f721fSFande Kong PetscInt ii, jj, ibeg, iend, irow; 1583c69f721fSFande Kong PetscInt *i, *j; 1584c69f721fSFande Kong PetscScalar *a; 1585c69f721fSFande Kong 1586c69f721fSFande Kong PetscFunctionBegin; 1587c69f721fSFande Kong if (!hA) PetscFunctionReturn(0); 1588c69f721fSFande Kong 158908defe43SFande Kong i = (PetscInt*) hypre_CSRMatrixI(hA); 159008defe43SFande Kong j = (PetscInt*) hypre_CSRMatrixJ(hA); 1591c69f721fSFande Kong a = hypre_CSRMatrixData(hA); 1592c69f721fSFande Kong 1593c69f721fSFande Kong for (ii = 0; ii < N; ii++) { 1594c69f721fSFande Kong irow = rows[ii]; 1595c69f721fSFande Kong ibeg = i[irow]; 1596c69f721fSFande Kong iend = i[irow+1]; 1597c69f721fSFande Kong for (jj = ibeg; jj < iend; jj++) 1598c69f721fSFande Kong if (j[jj] == irow) a[jj] = diag; 1599c69f721fSFande Kong else a[jj] = 0.0; 1600c69f721fSFande Kong } 1601c69f721fSFande Kong PetscFunctionReturn(0); 1602c69f721fSFande Kong } 1603c69f721fSFande Kong 1604ddbeb582SStefano Zampini static PetscErrorCode MatZeroRows_HYPRE(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b) 1605c69f721fSFande Kong { 1606c69f721fSFande Kong hypre_ParCSRMatrix *parcsr; 1607c69f721fSFande Kong PetscInt *lrows,len; 1608c69f721fSFande Kong PetscErrorCode ierr; 1609c69f721fSFande Kong 1610c69f721fSFande Kong PetscFunctionBegin; 1611c69f721fSFande Kong if (x || b) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Does not support to modify the solution and the right hand size \n"); 1612c69f721fSFande Kong /* retrieve the internal matrix */ 1613c69f721fSFande Kong ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr); 1614c69f721fSFande Kong /* get locally owned rows */ 1615c69f721fSFande Kong ierr = MatZeroRowsMapLocal_Private(A,N,rows,&len,&lrows);CHKERRQ(ierr); 1616c69f721fSFande Kong /* zero diagonal part */ 1617c69f721fSFande Kong ierr = MatZeroRows_HYPRE_CSRMatrix(hypre_ParCSRMatrixDiag(parcsr),len,lrows,diag);CHKERRQ(ierr); 1618c69f721fSFande Kong /* zero off-diagonal part */ 1619c69f721fSFande Kong ierr = MatZeroRows_HYPRE_CSRMatrix(hypre_ParCSRMatrixOffd(parcsr),len,lrows,0.0);CHKERRQ(ierr); 1620c69f721fSFande Kong 1621c69f721fSFande Kong ierr = PetscFree(lrows);CHKERRQ(ierr); 1622c69f721fSFande Kong PetscFunctionReturn(0); 1623c69f721fSFande Kong } 1624c69f721fSFande Kong 1625ddbeb582SStefano Zampini static PetscErrorCode MatAssemblyBegin_HYPRE(Mat mat,MatAssemblyType mode) 1626c69f721fSFande Kong { 1627c69f721fSFande Kong PetscErrorCode ierr; 1628c69f721fSFande Kong 1629c69f721fSFande Kong PetscFunctionBegin; 1630c69f721fSFande Kong if (mat->nooffprocentries) PetscFunctionReturn(0); 1631c69f721fSFande Kong 1632c69f721fSFande Kong ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr); 1633c69f721fSFande Kong PetscFunctionReturn(0); 1634c69f721fSFande Kong } 1635c69f721fSFande Kong 1636ddbeb582SStefano Zampini static PetscErrorCode MatGetRow_HYPRE(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 1637c69f721fSFande Kong { 1638c69f721fSFande Kong hypre_ParCSRMatrix *parcsr; 1639c69f721fSFande Kong PetscErrorCode ierr; 1640c69f721fSFande Kong 1641c69f721fSFande Kong PetscFunctionBegin; 1642c69f721fSFande Kong /* retrieve the internal matrix */ 1643c69f721fSFande Kong ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr); 1644c69f721fSFande Kong /* call HYPRE API */ 164508defe43SFande Kong PetscStackCallStandard(HYPRE_ParCSRMatrixGetRow,(parcsr,row,(HYPRE_Int*)nz,(HYPRE_Int**)idx,v)); 1646c69f721fSFande Kong PetscFunctionReturn(0); 1647c69f721fSFande Kong } 1648c69f721fSFande Kong 1649ddbeb582SStefano Zampini static PetscErrorCode MatRestoreRow_HYPRE(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 1650c69f721fSFande Kong { 1651c69f721fSFande Kong hypre_ParCSRMatrix *parcsr; 1652c69f721fSFande Kong PetscErrorCode ierr; 1653c69f721fSFande Kong 1654c69f721fSFande Kong PetscFunctionBegin; 1655c69f721fSFande Kong /* retrieve the internal matrix */ 1656c69f721fSFande Kong ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr); 1657c69f721fSFande Kong /* call HYPRE API */ 165808defe43SFande Kong PetscStackCallStandard(HYPRE_ParCSRMatrixRestoreRow,(parcsr,row,(HYPRE_Int*)nz,(HYPRE_Int**)idx,v)); 1659c69f721fSFande Kong PetscFunctionReturn(0); 1660c69f721fSFande Kong } 1661c69f721fSFande Kong 1662ddbeb582SStefano Zampini static PetscErrorCode MatGetValues_HYPRE(Mat A,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[]) 1663c69f721fSFande Kong { 166445b8d346SStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 1665c69f721fSFande Kong PetscInt i; 16661d4906efSStefano Zampini 1667c69f721fSFande Kong PetscFunctionBegin; 1668c69f721fSFande Kong if (!m || !n) PetscFunctionReturn(0); 1669c69f721fSFande Kong /* Ignore negative row indices 1670c69f721fSFande Kong * And negative column indices should be automatically ignored in hypre 1671c69f721fSFande Kong * */ 1672c69f721fSFande Kong for (i=0; i<m; i++) 167345b8d346SStefano Zampini if (idxm[i] >= 0) PetscStackCallStandard(HYPRE_IJMatrixGetValues,(hA->ij,1,(HYPRE_Int*)&n,(HYPRE_Int*)&idxm[i],(HYPRE_Int*)idxn,&v[i*n])); 1674c69f721fSFande Kong PetscFunctionReturn(0); 1675c69f721fSFande Kong } 1676c69f721fSFande Kong 1677ddbeb582SStefano Zampini static PetscErrorCode MatSetOption_HYPRE(Mat A,MatOption op,PetscBool flg) 1678ddbeb582SStefano Zampini { 1679ddbeb582SStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 1680ddbeb582SStefano Zampini 1681ddbeb582SStefano Zampini PetscFunctionBegin; 1682ddbeb582SStefano Zampini switch (op) 1683ddbeb582SStefano Zampini { 1684ddbeb582SStefano Zampini case MAT_NO_OFF_PROC_ENTRIES: 1685ddbeb582SStefano Zampini if (flg) { 1686ddbeb582SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixSetMaxOffProcElmts,(hA->ij,0)); 1687ddbeb582SStefano Zampini } 1688ddbeb582SStefano Zampini break; 1689ddbeb582SStefano Zampini default: 1690ddbeb582SStefano Zampini break; 1691ddbeb582SStefano Zampini } 1692ddbeb582SStefano Zampini PetscFunctionReturn(0); 1693ddbeb582SStefano Zampini } 1694c69f721fSFande Kong 169545b8d346SStefano Zampini static PetscErrorCode MatView_HYPRE(Mat A, PetscViewer view) 169645b8d346SStefano Zampini { 169745b8d346SStefano Zampini hypre_ParCSRMatrix *parcsr; 169845b8d346SStefano Zampini PetscErrorCode ierr; 169945b8d346SStefano Zampini Mat B; 170045b8d346SStefano Zampini PetscViewerFormat format; 170145b8d346SStefano Zampini PetscErrorCode (*mview)(Mat,PetscViewer) = NULL; 170245b8d346SStefano Zampini 170345b8d346SStefano Zampini PetscFunctionBegin; 170445b8d346SStefano Zampini ierr = PetscViewerGetFormat(view,&format);CHKERRQ(ierr); 170545b8d346SStefano Zampini if (format != PETSC_VIEWER_NATIVE) { 170645b8d346SStefano Zampini ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr); 170745b8d346SStefano Zampini ierr = MatCreateFromParCSR(parcsr,MATAIJ,PETSC_USE_POINTER,&B);CHKERRQ(ierr); 170845b8d346SStefano Zampini ierr = MatGetOperation(B,MATOP_VIEW,(void(**)(void))&mview);CHKERRQ(ierr); 170945b8d346SStefano Zampini if (!mview) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing view operation");CHKERRQ(ierr); 171045b8d346SStefano Zampini ierr = (*mview)(B,view);CHKERRQ(ierr); 171145b8d346SStefano Zampini ierr = MatDestroy(&B);CHKERRQ(ierr); 171245b8d346SStefano Zampini } else { 171345b8d346SStefano Zampini Mat_HYPRE *hA = (Mat_HYPRE*)A->data; 171445b8d346SStefano Zampini PetscMPIInt size; 171545b8d346SStefano Zampini PetscBool isascii; 171645b8d346SStefano Zampini const char *filename; 171745b8d346SStefano Zampini 171845b8d346SStefano Zampini /* HYPRE uses only text files */ 171945b8d346SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)view,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 172045b8d346SStefano Zampini if (!isascii) SETERRQ1(PetscObjectComm((PetscObject)view),PETSC_ERR_SUP,"PetscViewerType %s: native HYPRE format needs PETSCVIEWERASCII",((PetscObject)view)->type_name); 172145b8d346SStefano Zampini ierr = PetscViewerFileGetName(view,&filename);CHKERRQ(ierr); 172245b8d346SStefano Zampini PetscStackCallStandard(HYPRE_IJMatrixPrint,(hA->ij,filename)); 172345b8d346SStefano Zampini ierr = MPI_Comm_size(hA->comm,&size);CHKERRQ(ierr); 172445b8d346SStefano Zampini if (size > 1) { 172545b8d346SStefano Zampini ierr = PetscViewerASCIIPrintf(view,"Matrix files: %s.%05d ... %s.%05d\n",filename,0,filename,size-1);CHKERRQ(ierr); 172645b8d346SStefano Zampini } else { 172745b8d346SStefano Zampini ierr = PetscViewerASCIIPrintf(view,"Matrix file: %s.%05d\n",filename,0);CHKERRQ(ierr); 172845b8d346SStefano Zampini } 172945b8d346SStefano Zampini } 173045b8d346SStefano Zampini PetscFunctionReturn(0); 173145b8d346SStefano Zampini } 173245b8d346SStefano Zampini 173345b8d346SStefano Zampini static PetscErrorCode MatDuplicate_HYPRE(Mat A,MatDuplicateOption op, Mat *B) 173445b8d346SStefano Zampini { 173545b8d346SStefano Zampini hypre_ParCSRMatrix *parcsr; 173645b8d346SStefano Zampini PetscErrorCode ierr; 173745b8d346SStefano Zampini PetscCopyMode cpmode; 173845b8d346SStefano Zampini 173945b8d346SStefano Zampini PetscFunctionBegin; 174045b8d346SStefano Zampini ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr); 174145b8d346SStefano Zampini if (op == MAT_DO_NOT_COPY_VALUES || op == MAT_SHARE_NONZERO_PATTERN) { 174245b8d346SStefano Zampini parcsr = hypre_ParCSRMatrixCompleteClone(parcsr); 174345b8d346SStefano Zampini cpmode = PETSC_OWN_POINTER; 174445b8d346SStefano Zampini } else { 174545b8d346SStefano Zampini cpmode = PETSC_COPY_VALUES; 174645b8d346SStefano Zampini } 174745b8d346SStefano Zampini ierr = MatCreateFromParCSR(parcsr,MATHYPRE,cpmode,B);CHKERRQ(ierr); 174845b8d346SStefano Zampini PetscFunctionReturn(0); 174945b8d346SStefano Zampini } 175045b8d346SStefano Zampini 1751465edc17SStefano Zampini static PetscErrorCode MatCopy_HYPRE(Mat A, Mat B, MatStructure str) 1752465edc17SStefano Zampini { 1753465edc17SStefano Zampini hypre_ParCSRMatrix *acsr,*bcsr; 1754465edc17SStefano Zampini PetscErrorCode ierr; 1755465edc17SStefano Zampini 1756465edc17SStefano Zampini PetscFunctionBegin; 1757465edc17SStefano Zampini if (str == SAME_NONZERO_PATTERN && A->ops->copy == B->ops->copy) { 1758465edc17SStefano Zampini ierr = MatHYPREGetParCSR_HYPRE(A,&acsr);CHKERRQ(ierr); 1759465edc17SStefano Zampini ierr = MatHYPREGetParCSR_HYPRE(B,&bcsr);CHKERRQ(ierr); 1760465edc17SStefano Zampini PetscStackCallStandard(hypre_ParCSRMatrixCopy,(acsr,bcsr,1)); 1761465edc17SStefano Zampini ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1762465edc17SStefano Zampini ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1763465edc17SStefano Zampini } else { 1764465edc17SStefano Zampini ierr = MatCopy_Basic(A,B,str);CHKERRQ(ierr); 1765465edc17SStefano Zampini } 1766465edc17SStefano Zampini PetscFunctionReturn(0); 1767465edc17SStefano Zampini } 1768465edc17SStefano Zampini 17696305df00SStefano Zampini static PetscErrorCode MatGetDiagonal_HYPRE(Mat A, Vec d) 17706305df00SStefano Zampini { 17716305df00SStefano Zampini hypre_ParCSRMatrix *parcsr; 17726305df00SStefano Zampini hypre_CSRMatrix *dmat; 17736305df00SStefano Zampini PetscScalar *a,*data = NULL; 17746305df00SStefano Zampini PetscInt i,*diag = NULL; 17756305df00SStefano Zampini PetscBool cong; 17766305df00SStefano Zampini PetscErrorCode ierr; 17776305df00SStefano Zampini 17786305df00SStefano Zampini PetscFunctionBegin; 17796305df00SStefano Zampini ierr = MatHasCongruentLayouts(A,&cong);CHKERRQ(ierr); 17806305df00SStefano Zampini if (!cong) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only for square matrices with same local distributions of rows and columns"); 17816305df00SStefano Zampini #if defined(PETSC_USE_DEBUG) 17826305df00SStefano Zampini { 17836305df00SStefano Zampini PetscBool miss; 17846305df00SStefano Zampini ierr = MatMissingDiagonal(A,&miss,NULL);CHKERRQ(ierr); 17856305df00SStefano Zampini if (miss && A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented when diagonal entries are missing"); 17866305df00SStefano Zampini } 17876305df00SStefano Zampini #endif 17886305df00SStefano Zampini ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr); 17896305df00SStefano Zampini dmat = hypre_ParCSRMatrixDiag(parcsr); 17906305df00SStefano Zampini if (dmat) { 17916305df00SStefano Zampini ierr = VecGetArray(d,&a);CHKERRQ(ierr); 17926305df00SStefano Zampini diag = (PetscInt*)(hypre_CSRMatrixI(dmat)); 17936305df00SStefano Zampini data = (PetscScalar*)(hypre_CSRMatrixData(dmat)); 17946305df00SStefano Zampini for (i=0;i<A->rmap->n;i++) a[i] = data[diag[i]]; 17956305df00SStefano Zampini ierr = VecRestoreArray(d,&a);CHKERRQ(ierr); 17966305df00SStefano Zampini } 17976305df00SStefano Zampini PetscFunctionReturn(0); 17986305df00SStefano Zampini } 17996305df00SStefano Zampini 1800*363d496dSStefano Zampini #include <petscblaslapack.h> 1801*363d496dSStefano Zampini 1802*363d496dSStefano Zampini static PetscErrorCode MatAXPY_HYPRE(Mat Y,PetscScalar a,Mat X,MatStructure str) 1803*363d496dSStefano Zampini { 1804*363d496dSStefano Zampini PetscErrorCode ierr; 1805*363d496dSStefano Zampini 1806*363d496dSStefano Zampini PetscFunctionBegin; 1807*363d496dSStefano Zampini if (str == SAME_NONZERO_PATTERN) { 1808*363d496dSStefano Zampini hypre_ParCSRMatrix *x,*y; 1809*363d496dSStefano Zampini hypre_CSRMatrix *xloc,*yloc; 1810*363d496dSStefano Zampini PetscInt xnnz,ynnz; 1811*363d496dSStefano Zampini PetscScalar *xarr,*yarr; 1812*363d496dSStefano Zampini PetscBLASInt one=1,bnz; 1813*363d496dSStefano Zampini 1814*363d496dSStefano Zampini ierr = MatHYPREGetParCSR(Y,&y);CHKERRQ(ierr); 1815*363d496dSStefano Zampini ierr = MatHYPREGetParCSR(X,&x);CHKERRQ(ierr); 1816*363d496dSStefano Zampini 1817*363d496dSStefano Zampini /* diagonal block */ 1818*363d496dSStefano Zampini xloc = hypre_ParCSRMatrixDiag(x); 1819*363d496dSStefano Zampini yloc = hypre_ParCSRMatrixDiag(y); 1820*363d496dSStefano Zampini xnnz = 0; 1821*363d496dSStefano Zampini ynnz = 0; 1822*363d496dSStefano Zampini xarr = NULL; 1823*363d496dSStefano Zampini yarr = NULL; 1824*363d496dSStefano Zampini if (xloc) { 1825*363d496dSStefano Zampini xarr = (PetscScalar*)(hypre_CSRMatrixData(xloc)); 1826*363d496dSStefano Zampini xnnz = hypre_CSRMatrixNumNonzeros(xloc); 1827*363d496dSStefano Zampini } 1828*363d496dSStefano Zampini if (yloc) { 1829*363d496dSStefano Zampini yarr = (PetscScalar*)(hypre_CSRMatrixData(yloc)); 1830*363d496dSStefano Zampini ynnz = hypre_CSRMatrixNumNonzeros(yloc); 1831*363d496dSStefano Zampini } 1832*363d496dSStefano Zampini if (xnnz != ynnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Different number of nonzeros in diagonal block %D != %D",xnnz,ynnz); 1833*363d496dSStefano Zampini ierr = PetscBLASIntCast(xnnz,&bnz);CHKERRQ(ierr); 1834*363d496dSStefano Zampini PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&a,xarr,&one,yarr,&one)); 1835*363d496dSStefano Zampini 1836*363d496dSStefano Zampini /* off-diagonal block */ 1837*363d496dSStefano Zampini xloc = hypre_ParCSRMatrixOffd(x); 1838*363d496dSStefano Zampini yloc = hypre_ParCSRMatrixOffd(y); 1839*363d496dSStefano Zampini xnnz = 0; 1840*363d496dSStefano Zampini ynnz = 0; 1841*363d496dSStefano Zampini xarr = NULL; 1842*363d496dSStefano Zampini yarr = NULL; 1843*363d496dSStefano Zampini if (xloc) { 1844*363d496dSStefano Zampini xarr = (PetscScalar*)(hypre_CSRMatrixData(xloc)); 1845*363d496dSStefano Zampini xnnz = hypre_CSRMatrixNumNonzeros(xloc); 1846*363d496dSStefano Zampini } 1847*363d496dSStefano Zampini if (yloc) { 1848*363d496dSStefano Zampini yarr = (PetscScalar*)(hypre_CSRMatrixData(yloc)); 1849*363d496dSStefano Zampini ynnz = hypre_CSRMatrixNumNonzeros(yloc); 1850*363d496dSStefano Zampini } 1851*363d496dSStefano Zampini if (xnnz != ynnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Different number of nonzeros in off-diagonal block %D != %D",xnnz,ynnz); 1852*363d496dSStefano Zampini ierr = PetscBLASIntCast(xnnz,&bnz);CHKERRQ(ierr); 1853*363d496dSStefano Zampini PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&bnz,&a,xarr,&one,yarr,&one)); 1854*363d496dSStefano Zampini } else if (str == SUBSET_NONZERO_PATTERN) { 1855*363d496dSStefano Zampini ierr = MatAXPY_Basic(Y,a,X,str);CHKERRQ(ierr); 1856*363d496dSStefano Zampini } else { 1857*363d496dSStefano Zampini Mat B; 1858*363d496dSStefano Zampini 1859*363d496dSStefano Zampini ierr = MatAXPY_Basic_Preallocate(Y,X,&B);CHKERRQ(ierr); 1860*363d496dSStefano Zampini ierr = MatAXPY_BasicWithPreallocation(B,Y,a,X,str);CHKERRQ(ierr); 1861*363d496dSStefano Zampini ierr = MatHeaderReplace(Y,&B);CHKERRQ(ierr); 1862*363d496dSStefano Zampini } 1863*363d496dSStefano Zampini PetscFunctionReturn(0); 1864*363d496dSStefano Zampini } 1865*363d496dSStefano Zampini 1866a055b5aaSBarry Smith /*MC 1867a055b5aaSBarry Smith MATHYPRE - MATHYPRE = "hypre" - A matrix type to be used for sequential and parallel sparse matrices 1868a055b5aaSBarry Smith based on the hypre IJ interface. 1869a055b5aaSBarry Smith 1870a055b5aaSBarry Smith Level: intermediate 1871a055b5aaSBarry Smith 1872a055b5aaSBarry Smith .seealso: MatCreate() 1873a055b5aaSBarry Smith M*/ 1874a055b5aaSBarry Smith 187563c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B) 187663c07aadSStefano Zampini { 187763c07aadSStefano Zampini Mat_HYPRE *hB; 187863c07aadSStefano Zampini PetscErrorCode ierr; 187963c07aadSStefano Zampini 188063c07aadSStefano Zampini PetscFunctionBegin; 188163c07aadSStefano Zampini ierr = PetscNewLog(B,&hB);CHKERRQ(ierr); 1882978814f1SStefano Zampini hB->inner_free = PETSC_TRUE; 1883c69f721fSFande Kong hB->available = PETSC_TRUE; 1884c69f721fSFande Kong hB->size = 0; 1885c69f721fSFande Kong hB->array = NULL; 1886978814f1SStefano Zampini 188763c07aadSStefano Zampini B->data = (void*)hB; 188863c07aadSStefano Zampini B->rmap->bs = 1; 188963c07aadSStefano Zampini B->assembled = PETSC_FALSE; 189063c07aadSStefano Zampini 1891de393100SStefano Zampini ierr = PetscMemzero(B->ops,sizeof(struct _MatOps));CHKERRQ(ierr); 189263c07aadSStefano Zampini B->ops->mult = MatMult_HYPRE; 189363c07aadSStefano Zampini B->ops->multtranspose = MatMultTranspose_HYPRE; 1894414bd5c3SStefano Zampini B->ops->multadd = MatMultAdd_HYPRE; 1895414bd5c3SStefano Zampini B->ops->multtransposeadd = MatMultTransposeAdd_HYPRE; 189663c07aadSStefano Zampini B->ops->setup = MatSetUp_HYPRE; 189763c07aadSStefano Zampini B->ops->destroy = MatDestroy_HYPRE; 189863c07aadSStefano Zampini B->ops->assemblyend = MatAssemblyEnd_HYPRE; 1899c69f721fSFande Kong B->ops->assemblybegin = MatAssemblyBegin_HYPRE; 1900cd8bc7baSStefano Zampini B->ops->ptap = MatPtAP_HYPRE_HYPRE; 1901d501dc42Sstefano_zampini B->ops->matmult = MatMatMult_HYPRE_HYPRE; 1902d975228cSstefano_zampini B->ops->setvalues = MatSetValues_HYPRE; 190368ec7858SStefano Zampini B->ops->missingdiagonal = MatMissingDiagonal_HYPRE; 190468ec7858SStefano Zampini B->ops->scale = MatScale_HYPRE; 190568ec7858SStefano Zampini B->ops->zerorowscolumns = MatZeroRowsColumns_HYPRE; 1906c69f721fSFande Kong B->ops->zeroentries = MatZeroEntries_HYPRE; 1907c69f721fSFande Kong B->ops->zerorows = MatZeroRows_HYPRE; 1908c69f721fSFande Kong B->ops->getrow = MatGetRow_HYPRE; 1909c69f721fSFande Kong B->ops->restorerow = MatRestoreRow_HYPRE; 1910c69f721fSFande Kong B->ops->getvalues = MatGetValues_HYPRE; 1911ddbeb582SStefano Zampini B->ops->setoption = MatSetOption_HYPRE; 191245b8d346SStefano Zampini B->ops->duplicate = MatDuplicate_HYPRE; 1913465edc17SStefano Zampini B->ops->copy = MatCopy_HYPRE; 191445b8d346SStefano Zampini B->ops->view = MatView_HYPRE; 19156305df00SStefano Zampini B->ops->getdiagonal = MatGetDiagonal_HYPRE; 1916*363d496dSStefano Zampini B->ops->axpy = MatAXPY_HYPRE; 191745b8d346SStefano Zampini 191845b8d346SStefano Zampini /* build cache for off array entries formed */ 191945b8d346SStefano Zampini ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);CHKERRQ(ierr); 192063c07aadSStefano Zampini 192163c07aadSStefano Zampini ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr); 192263c07aadSStefano Zampini ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr); 192363c07aadSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr); 19242df22349SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr); 1925c1a070e6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr); 1926c1a070e6SStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr); 1927d975228cSstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPRESetPreallocation_C",MatHYPRESetPreallocation_HYPRE);CHKERRQ(ierr); 1928dd9c0a25Sstefano_zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPREGetParCSR_C",MatHYPREGetParCSR_HYPRE);CHKERRQ(ierr); 192975d48cdbSStefano Zampini ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_is_hypre_C",MatPtAP_IS_XAIJ);CHKERRQ(ierr); 193063c07aadSStefano Zampini PetscFunctionReturn(0); 193163c07aadSStefano Zampini } 193263c07aadSStefano Zampini 1933225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void *ptr) 1934225daaf8SStefano Zampini { 1935225daaf8SStefano Zampini PetscFunctionBegin; 1936e6de0934SSatish Balay hypre_TFree(ptr,HYPRE_MEMORY_HOST); 1937225daaf8SStefano Zampini PetscFunctionReturn(0); 1938225daaf8SStefano Zampini } 1939