xref: /petsc/src/mat/impls/hypre/mhypre.c (revision dd9c0a25dbe2a535db1b669ab1f65e7c9e227acd)
163c07aadSStefano Zampini 
263c07aadSStefano Zampini /*
363c07aadSStefano Zampini     Creates hypre ijmatrix from PETSc matrix
463c07aadSStefano Zampini */
5225daaf8SStefano Zampini 
6225daaf8SStefano Zampini /*MC
7225daaf8SStefano Zampini    MATHYPRE - MATHYPRE = "hypre" - A matrix type to be used for sequential and parallel sparse matrices
8225daaf8SStefano Zampini           based on the hypre IJ interface.
9225daaf8SStefano Zampini 
10225daaf8SStefano Zampini    Level: intermediate
11225daaf8SStefano Zampini 
12225daaf8SStefano Zampini .seealso: MatCreate()
13225daaf8SStefano Zampini M*/
14225daaf8SStefano Zampini 
15*dd9c0a25Sstefano_zampini #include <petscmathypre.h>
1663c07aadSStefano Zampini #include <petsc/private/matimpl.h>
1763c07aadSStefano Zampini #include <../src/mat/impls/hypre/mhypre.h>
1863c07aadSStefano Zampini #include <../src/mat/impls/aij/mpi/mpiaij.h>
1958968eb6SStefano Zampini #include <../src/vec/vec/impls/hypre/vhyp.h>
2058968eb6SStefano Zampini #include <HYPRE.h>
21c1a070e6SStefano Zampini #include <HYPRE_utilities.h>
22cd8bc7baSStefano Zampini #include <_hypre_parcsr_ls.h>
2363c07aadSStefano Zampini 
2463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*);
2563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix);
2663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
2763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
2863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool);
29225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void*);
3063c07aadSStefano Zampini 
3163c07aadSStefano Zampini #undef __FUNCT__
3263c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate"
3363c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij)
3463c07aadSStefano Zampini {
3563c07aadSStefano Zampini   PetscErrorCode ierr;
3663c07aadSStefano Zampini   PetscInt       i,n_d,n_o;
3763c07aadSStefano Zampini   const PetscInt *ia_d,*ia_o;
3863c07aadSStefano Zampini   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
3963c07aadSStefano Zampini   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
4063c07aadSStefano Zampini 
4163c07aadSStefano Zampini   PetscFunctionBegin;
4263c07aadSStefano Zampini   if (A_d) { /* determine number of nonzero entries in local diagonal part */
4363c07aadSStefano Zampini     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
4463c07aadSStefano Zampini     if (done_d) {
4563c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr);
4663c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
4763c07aadSStefano Zampini         nnz_d[i] = ia_d[i+1] - ia_d[i];
4863c07aadSStefano Zampini       }
4963c07aadSStefano Zampini     }
5063c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
5163c07aadSStefano Zampini   }
5263c07aadSStefano Zampini   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
5363c07aadSStefano Zampini     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
5463c07aadSStefano Zampini     if (done_o) {
5563c07aadSStefano Zampini       ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr);
5663c07aadSStefano Zampini       for (i=0; i<n_o; i++) {
5763c07aadSStefano Zampini         nnz_o[i] = ia_o[i+1] - ia_o[i];
5863c07aadSStefano Zampini       }
5963c07aadSStefano Zampini     }
6063c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
6163c07aadSStefano Zampini   }
6263c07aadSStefano Zampini   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
6363c07aadSStefano Zampini     if (!done_o) { /* only diagonal part */
6463c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
6563c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
6663c07aadSStefano Zampini         nnz_o[i] = 0;
6763c07aadSStefano Zampini       }
6863c07aadSStefano Zampini     }
6963c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o));
7063c07aadSStefano Zampini     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
7163c07aadSStefano Zampini     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
7263c07aadSStefano Zampini   }
7363c07aadSStefano Zampini   PetscFunctionReturn(0);
7463c07aadSStefano Zampini }
7563c07aadSStefano Zampini 
7663c07aadSStefano Zampini #undef __FUNCT__
7763c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_CreateFromMat"
7863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA)
7963c07aadSStefano Zampini {
8063c07aadSStefano Zampini   PetscErrorCode ierr;
8163c07aadSStefano Zampini   PetscInt       rstart,rend,cstart,cend;
8263c07aadSStefano Zampini 
8363c07aadSStefano Zampini   PetscFunctionBegin;
84d975228cSstefano_zampini   ierr   = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
85d975228cSstefano_zampini   ierr   = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
8663c07aadSStefano Zampini   rstart = A->rmap->rstart;
8763c07aadSStefano Zampini   rend   = A->rmap->rend;
8863c07aadSStefano Zampini   cstart = A->cmap->rstart;
8963c07aadSStefano Zampini   cend   = A->cmap->rend;
9063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
9163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
9263c07aadSStefano Zampini   {
9363c07aadSStefano Zampini     PetscBool      same;
9463c07aadSStefano Zampini     Mat            A_d,A_o;
9563c07aadSStefano Zampini     const PetscInt *colmap;
9663c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
9763c07aadSStefano Zampini     if (same) {
9863c07aadSStefano Zampini       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
9963c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
10063c07aadSStefano Zampini       PetscFunctionReturn(0);
10163c07aadSStefano Zampini     }
10263c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
10363c07aadSStefano Zampini     if (same) {
10463c07aadSStefano Zampini       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
10563c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
10663c07aadSStefano Zampini       PetscFunctionReturn(0);
10763c07aadSStefano Zampini     }
10863c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
10963c07aadSStefano Zampini     if (same) {
11063c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
11163c07aadSStefano Zampini       PetscFunctionReturn(0);
11263c07aadSStefano Zampini     }
11363c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
11463c07aadSStefano Zampini     if (same) {
11563c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
11663c07aadSStefano Zampini       PetscFunctionReturn(0);
11763c07aadSStefano Zampini     }
11863c07aadSStefano Zampini   }
11963c07aadSStefano Zampini   PetscFunctionReturn(0);
12063c07aadSStefano Zampini }
12163c07aadSStefano Zampini 
12263c07aadSStefano Zampini #undef __FUNCT__
12363c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixCopy"
12463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij)
12563c07aadSStefano Zampini {
12663c07aadSStefano Zampini   PetscErrorCode    ierr;
12763c07aadSStefano Zampini   PetscInt          i,rstart,rend,ncols,nr,nc;
12863c07aadSStefano Zampini   const PetscScalar *values;
12963c07aadSStefano Zampini   const PetscInt    *cols;
13063c07aadSStefano Zampini   PetscBool         flg;
13163c07aadSStefano Zampini 
13263c07aadSStefano Zampini   PetscFunctionBegin;
13363c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
13463c07aadSStefano Zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
13563c07aadSStefano Zampini   if (flg && nr == nc) {
13663c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
13763c07aadSStefano Zampini     PetscFunctionReturn(0);
13863c07aadSStefano Zampini   }
13963c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
14063c07aadSStefano Zampini   if (flg) {
14163c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
14263c07aadSStefano Zampini     PetscFunctionReturn(0);
14363c07aadSStefano Zampini   }
14463c07aadSStefano Zampini 
14563c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
14663c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
14763c07aadSStefano Zampini   for (i=rstart; i<rend; i++) {
14863c07aadSStefano Zampini     ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
14963c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
15063c07aadSStefano Zampini     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
15163c07aadSStefano Zampini   }
15263c07aadSStefano Zampini   PetscFunctionReturn(0);
15363c07aadSStefano Zampini }
15463c07aadSStefano Zampini 
15563c07aadSStefano Zampini #undef __FUNCT__
15663c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ"
15763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij)
15863c07aadSStefano Zampini {
15963c07aadSStefano Zampini   PetscErrorCode        ierr;
16063c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag = (Mat_SeqAIJ*)A->data;
16158968eb6SStefano Zampini   HYPRE_Int             type;
16263c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
16363c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
16463c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag;
16563c07aadSStefano Zampini 
16663c07aadSStefano Zampini   PetscFunctionBegin;
16763c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
168ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
169ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
170ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
17163c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
17263c07aadSStefano Zampini   /*
17363c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
17463c07aadSStefano Zampini   */
17563c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));
17663c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));
17763c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
178ea9daf28SStefano Zampini 
179ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
18063c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
18163c07aadSStefano Zampini   PetscFunctionReturn(0);
18263c07aadSStefano Zampini }
18363c07aadSStefano Zampini 
18463c07aadSStefano Zampini #undef __FUNCT__
18563c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ"
18663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
18763c07aadSStefano Zampini {
18863c07aadSStefano Zampini   PetscErrorCode        ierr;
18963c07aadSStefano Zampini   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
19063c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag,*poffd;
19163c07aadSStefano Zampini   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
19258968eb6SStefano Zampini   HYPRE_Int             type;
19363c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
19463c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
19563c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag,*hoffd;
19663c07aadSStefano Zampini 
19763c07aadSStefano Zampini   PetscFunctionBegin;
19863c07aadSStefano Zampini   pdiag = (Mat_SeqAIJ*) pA->A->data;
19963c07aadSStefano Zampini   poffd = (Mat_SeqAIJ*) pA->B->data;
20063c07aadSStefano Zampini   /* cstart is only valid for square MPIAIJ layed out in the usual way */
20163c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
20263c07aadSStefano Zampini 
20363c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
204ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
205ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
206ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
20763c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
20863c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
20963c07aadSStefano Zampini 
21063c07aadSStefano Zampini   /*
21163c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
21263c07aadSStefano Zampini   */
21363c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
21463c07aadSStefano Zampini   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
21563c07aadSStefano Zampini   jj  = (PetscInt*)hdiag->j;
21663c07aadSStefano Zampini   pjj = (PetscInt*)pdiag->j;
21763c07aadSStefano Zampini   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
21863c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
21963c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
22063c07aadSStefano Zampini   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
22163c07aadSStefano Zampini      If we hacked a hypre a bit more we might be able to avoid this step */
22263c07aadSStefano Zampini   jj  = (PetscInt*) hoffd->j;
22363c07aadSStefano Zampini   pjj = (PetscInt*) poffd->j;
22463c07aadSStefano Zampini   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
22563c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
22663c07aadSStefano Zampini 
227ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
22863c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
22963c07aadSStefano Zampini   PetscFunctionReturn(0);
23063c07aadSStefano Zampini }
23163c07aadSStefano Zampini 
23263c07aadSStefano Zampini #undef __FUNCT__
2332df22349SStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_IS"
2342df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B)
2352df22349SStefano Zampini {
2362df22349SStefano Zampini   Mat_HYPRE*             mhA = (Mat_HYPRE*)(A->data);
2372df22349SStefano Zampini   Mat                    lA;
2382df22349SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
2392df22349SStefano Zampini   IS                     is;
2402df22349SStefano Zampini   hypre_ParCSRMatrix     *hA;
2412df22349SStefano Zampini   hypre_CSRMatrix        *hdiag,*hoffd;
2422df22349SStefano Zampini   MPI_Comm               comm;
2432df22349SStefano Zampini   PetscScalar            *hdd,*hod,*aa,*data;
2442df22349SStefano Zampini   PetscInt               *col_map_offd,*hdi,*hdj,*hoi,*hoj;
2452df22349SStefano Zampini   PetscInt               *ii,*jj,*iptr,*jptr;
2462df22349SStefano Zampini   PetscInt               cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N;
24758968eb6SStefano Zampini   HYPRE_Int              type;
2482df22349SStefano Zampini   PetscErrorCode         ierr;
2492df22349SStefano Zampini 
2502df22349SStefano Zampini   PetscFunctionBegin;
251a1787963SStefano Zampini   comm = PetscObjectComm((PetscObject)A);
2522df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type));
2532df22349SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
2542df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA));
2552df22349SStefano Zampini   M     = hypre_ParCSRMatrixGlobalNumRows(hA);
2562df22349SStefano Zampini   N     = hypre_ParCSRMatrixGlobalNumCols(hA);
2572df22349SStefano Zampini   str   = hypre_ParCSRMatrixFirstRowIndex(hA);
2582df22349SStefano Zampini   stc   = hypre_ParCSRMatrixFirstColDiag(hA);
2592df22349SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(hA);
2602df22349SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(hA);
2612df22349SStefano Zampini   dr    = hypre_CSRMatrixNumRows(hdiag);
2622df22349SStefano Zampini   dc    = hypre_CSRMatrixNumCols(hdiag);
2632df22349SStefano Zampini   nnz   = hypre_CSRMatrixNumNonzeros(hdiag);
2642df22349SStefano Zampini   hdi   = hypre_CSRMatrixI(hdiag);
2652df22349SStefano Zampini   hdj   = hypre_CSRMatrixJ(hdiag);
2662df22349SStefano Zampini   hdd   = hypre_CSRMatrixData(hdiag);
2672df22349SStefano Zampini   oc    = hypre_CSRMatrixNumCols(hoffd);
2682df22349SStefano Zampini   nnz  += hypre_CSRMatrixNumNonzeros(hoffd);
2692df22349SStefano Zampini   hoi   = hypre_CSRMatrixI(hoffd);
2702df22349SStefano Zampini   hoj   = hypre_CSRMatrixJ(hoffd);
2712df22349SStefano Zampini   hod   = hypre_CSRMatrixData(hoffd);
2722df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
2732df22349SStefano Zampini     PetscInt *aux;
2742df22349SStefano Zampini 
2752df22349SStefano Zampini     /* generate l2g maps for rows and cols */
2762df22349SStefano Zampini     ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr);
2772df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
2782df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2792df22349SStefano Zampini     col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
2802df22349SStefano Zampini     ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr);
2812df22349SStefano Zampini     for (i=0; i<dc; i++) aux[i] = i+stc;
2822df22349SStefano Zampini     for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i];
2832df22349SStefano Zampini     ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
2842df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
2852df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2862df22349SStefano Zampini     /* create MATIS object */
2872df22349SStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
2882df22349SStefano Zampini     ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr);
2892df22349SStefano Zampini     ierr = MatSetType(*B,MATIS);CHKERRQ(ierr);
2902df22349SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr);
2912df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
2922df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
2932df22349SStefano Zampini 
2942df22349SStefano Zampini     /* allocate CSR for local matrix */
2952df22349SStefano Zampini     ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr);
2962df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr);
2972df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr);
2982df22349SStefano Zampini   } else {
2992df22349SStefano Zampini     PetscInt  nr;
3002df22349SStefano Zampini     PetscBool done;
3012df22349SStefano Zampini     ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr);
3022df22349SStefano Zampini     ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr);
3032df22349SStefano Zampini     if (nr != dr) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of rows in local mat! %D != %D",nr,dr);
3042df22349SStefano Zampini     if (iptr[nr] < nnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros in local mat! reuse %D requested %D",iptr[nr],nnz);
3052df22349SStefano Zampini     ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr);
3062df22349SStefano Zampini   }
3072df22349SStefano Zampini   /* merge local matrices */
3082df22349SStefano Zampini   ii   = iptr;
3092df22349SStefano Zampini   jj   = jptr;
3102df22349SStefano Zampini   aa   = data;
3112df22349SStefano Zampini   *ii  = *(hdi++) + *(hoi++);
3122df22349SStefano Zampini   for (jd=0,jo=0,cum=0; *ii<nnz; cum++) {
3132df22349SStefano Zampini     PetscScalar *aold = aa;
3142df22349SStefano Zampini     PetscInt    *jold = jj,nc = jd+jo;
3152df22349SStefano Zampini     for (; jd<*hdi; jd++) { *jj++ = *hdj++;      *aa++ = *hdd++; }
3162df22349SStefano Zampini     for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
3172df22349SStefano Zampini     *(++ii) = *(hdi++) + *(hoi++);
3182df22349SStefano Zampini     ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr);
3192df22349SStefano Zampini   }
3202df22349SStefano Zampini   for (; cum<dr; cum++) *(++ii) = nnz;
3212df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
322a033916dSStefano Zampini     Mat_SeqAIJ* a;
323a033916dSStefano Zampini 
3242df22349SStefano Zampini     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr);
3252df22349SStefano Zampini     ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr);
326a033916dSStefano Zampini     /* hack SeqAIJ */
327a033916dSStefano Zampini     a          = (Mat_SeqAIJ*)(lA->data);
328a033916dSStefano Zampini     a->free_a  = PETSC_TRUE;
329a033916dSStefano Zampini     a->free_ij = PETSC_TRUE;
3302df22349SStefano Zampini     ierr = MatDestroy(&lA);CHKERRQ(ierr);
3312df22349SStefano Zampini   }
3322df22349SStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3332df22349SStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3342df22349SStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
3352df22349SStefano Zampini     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
3362df22349SStefano Zampini   }
3372df22349SStefano Zampini   PetscFunctionReturn(0);
3382df22349SStefano Zampini }
3392df22349SStefano Zampini 
3402df22349SStefano Zampini #undef __FUNCT__
34163c07aadSStefano Zampini #define __FUNCT__ "MatConvert_AIJ_HYPRE"
34263c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
34363c07aadSStefano Zampini {
34463c07aadSStefano Zampini   Mat_HYPRE      *hB;
34563c07aadSStefano Zampini   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
34663c07aadSStefano Zampini   PetscErrorCode ierr;
34763c07aadSStefano Zampini 
34863c07aadSStefano Zampini   PetscFunctionBegin;
34963c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
35063c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
35163c07aadSStefano Zampini     /* always destroy the old matrix and create a new memory;
35263c07aadSStefano Zampini        hope this does not churn the memory too much. The problem
35363c07aadSStefano Zampini        is I do not know if it is possible to put the matrix back to
35463c07aadSStefano Zampini        its initial state so that we can directly copy the values
35563c07aadSStefano Zampini        the second time through. */
35663c07aadSStefano Zampini     hB = (Mat_HYPRE*)((*B)->data);
35763c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
35863c07aadSStefano Zampini   } else {
35963c07aadSStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
36063c07aadSStefano Zampini     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
36163c07aadSStefano Zampini     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
36263c07aadSStefano Zampini     ierr = MatSetUp(*B);CHKERRQ(ierr);
36363c07aadSStefano Zampini     hB   = (Mat_HYPRE*)((*B)->data);
36463c07aadSStefano Zampini   }
36563c07aadSStefano Zampini   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
36663c07aadSStefano Zampini   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
36763c07aadSStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
36863c07aadSStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
36963c07aadSStefano Zampini   PetscFunctionReturn(0);
37063c07aadSStefano Zampini }
37163c07aadSStefano Zampini 
37263c07aadSStefano Zampini #undef __FUNCT__
37363c07aadSStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_AIJ"
374ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
37563c07aadSStefano Zampini {
37663c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
37763c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
37863c07aadSStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
37963c07aadSStefano Zampini   MPI_Comm           comm;
38063c07aadSStefano Zampini   PetscScalar        *da,*oa,*aptr;
38163c07aadSStefano Zampini   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
38263c07aadSStefano Zampini   PetscInt           i,dnnz,onnz,m,n;
38358968eb6SStefano Zampini   HYPRE_Int          type;
38463c07aadSStefano Zampini   PetscMPIInt        size;
38563c07aadSStefano Zampini   PetscErrorCode     ierr;
38663c07aadSStefano Zampini 
38763c07aadSStefano Zampini   PetscFunctionBegin;
38863c07aadSStefano Zampini   comm = PetscObjectComm((PetscObject)A);
38963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
39063c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
39163c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
39263c07aadSStefano Zampini     PetscBool ismpiaij,isseqaij;
39363c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
39463c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
39563c07aadSStefano Zampini     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
39663c07aadSStefano Zampini   }
39763c07aadSStefano Zampini   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
39863c07aadSStefano Zampini 
39963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
40063c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(parcsr);
40163c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(parcsr);
40263c07aadSStefano Zampini   m     = hypre_CSRMatrixNumRows(hdiag);
40363c07aadSStefano Zampini   n     = hypre_CSRMatrixNumCols(hdiag);
40463c07aadSStefano Zampini   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
40563c07aadSStefano Zampini   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
406225daaf8SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
40763c07aadSStefano Zampini     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
40863c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
40963c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
410225daaf8SStefano Zampini   } else if (reuse == MAT_REUSE_MATRIX) {
41163c07aadSStefano Zampini     PetscInt  nr;
41263c07aadSStefano Zampini     PetscBool done;
41363c07aadSStefano Zampini     if (size > 1) {
41463c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
41563c07aadSStefano Zampini 
41663c07aadSStefano Zampini       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
41763c07aadSStefano Zampini       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows in diag part! %D != %D",nr,m);
41863c07aadSStefano Zampini       if (dii[nr] < dnnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros in diag part! reuse %D hypre %D",dii[nr],dnnz);
41963c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
42063c07aadSStefano Zampini     } else {
42163c07aadSStefano Zampini       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
42263c07aadSStefano Zampini       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
42363c07aadSStefano Zampini       if (dii[nr] < dnnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros! reuse %D hypre %D",dii[nr],dnnz);
42463c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
42563c07aadSStefano Zampini     }
426225daaf8SStefano Zampini   } else { /* MAT_INPLACE_MATRIX */
427225daaf8SStefano Zampini     dii = hypre_CSRMatrixI(hdiag);
428225daaf8SStefano Zampini     djj = hypre_CSRMatrixJ(hdiag);
429225daaf8SStefano Zampini     da  = hypre_CSRMatrixData(hdiag);
43063c07aadSStefano Zampini   }
43163c07aadSStefano Zampini   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
43263c07aadSStefano Zampini   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
43363c07aadSStefano Zampini   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
43463c07aadSStefano Zampini   iptr = djj;
43563c07aadSStefano Zampini   aptr = da;
43663c07aadSStefano Zampini   for (i=0; i<m; i++) {
43763c07aadSStefano Zampini     PetscInt nc = dii[i+1]-dii[i];
43863c07aadSStefano Zampini     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
43963c07aadSStefano Zampini     iptr += nc;
44063c07aadSStefano Zampini     aptr += nc;
44163c07aadSStefano Zampini   }
44263c07aadSStefano Zampini   if (size > 1) {
44363c07aadSStefano Zampini     PetscInt *offdj,*coffd;
44463c07aadSStefano Zampini 
445225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
44663c07aadSStefano Zampini       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
44763c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
44863c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
449225daaf8SStefano Zampini     } else if (reuse == MAT_REUSE_MATRIX) {
45063c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
45163c07aadSStefano Zampini       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
45263c07aadSStefano Zampini       PetscBool  done;
45363c07aadSStefano Zampini 
45463c07aadSStefano Zampini       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
45563c07aadSStefano Zampini       if (nr != hr) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows in offdiag part! %D != %D",nr,hr);
45663c07aadSStefano Zampini       if (oii[nr] < onnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros in offdiag part! reuse %D hypre %D",oii[nr],onnz);
45763c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
458225daaf8SStefano Zampini     } else { /* MAT_INPLACE_MATRIX */
459225daaf8SStefano Zampini       oii = hypre_CSRMatrixI(hoffd);
460225daaf8SStefano Zampini       ojj = hypre_CSRMatrixJ(hoffd);
461225daaf8SStefano Zampini       oa  = hypre_CSRMatrixData(hoffd);
46263c07aadSStefano Zampini     }
46363c07aadSStefano Zampini     ierr  = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
46463c07aadSStefano Zampini     offdj = hypre_CSRMatrixJ(hoffd);
46563c07aadSStefano Zampini     coffd = hypre_ParCSRMatrixColMapOffd(parcsr);
46663c07aadSStefano Zampini     for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]];
46763c07aadSStefano Zampini     ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr);
46863c07aadSStefano Zampini     iptr = ojj;
46963c07aadSStefano Zampini     aptr = oa;
47063c07aadSStefano Zampini     for (i=0; i<m; i++) {
47163c07aadSStefano Zampini        PetscInt nc = oii[i+1]-oii[i];
47263c07aadSStefano Zampini        ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
47363c07aadSStefano Zampini        iptr += nc;
47463c07aadSStefano Zampini        aptr += nc;
47563c07aadSStefano Zampini     }
476225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
47763c07aadSStefano Zampini       Mat_MPIAIJ *b;
47863c07aadSStefano Zampini       Mat_SeqAIJ *d,*o;
479225daaf8SStefano Zampini 
48041073d71Sstefano_zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,B);CHKERRQ(ierr);
48163c07aadSStefano Zampini       /* hack MPIAIJ */
48263c07aadSStefano Zampini       b          = (Mat_MPIAIJ*)((*B)->data);
48363c07aadSStefano Zampini       d          = (Mat_SeqAIJ*)b->A->data;
48463c07aadSStefano Zampini       o          = (Mat_SeqAIJ*)b->B->data;
48563c07aadSStefano Zampini       d->free_a  = PETSC_TRUE;
48663c07aadSStefano Zampini       d->free_ij = PETSC_TRUE;
48763c07aadSStefano Zampini       o->free_a  = PETSC_TRUE;
48863c07aadSStefano Zampini       o->free_ij = PETSC_TRUE;
489225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
490225daaf8SStefano Zampini       Mat T;
49141073d71Sstefano_zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,&T);CHKERRQ(ierr);
492225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
493225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
494225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
495225daaf8SStefano Zampini       hypre_CSRMatrixI(hoffd)    = NULL;
496225daaf8SStefano Zampini       hypre_CSRMatrixJ(hoffd)    = NULL;
497225daaf8SStefano Zampini       hypre_CSRMatrixData(hoffd) = NULL;
498225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
49963c07aadSStefano Zampini     }
500225daaf8SStefano Zampini   } else {
501225daaf8SStefano Zampini     oii  = NULL;
502225daaf8SStefano Zampini     ojj  = NULL;
503225daaf8SStefano Zampini     oa   = NULL;
504225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
50563c07aadSStefano Zampini       Mat_SeqAIJ* b;
50663c07aadSStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
50763c07aadSStefano Zampini       /* hack SeqAIJ */
50863c07aadSStefano Zampini       b          = (Mat_SeqAIJ*)((*B)->data);
50963c07aadSStefano Zampini       b->free_a  = PETSC_TRUE;
51063c07aadSStefano Zampini       b->free_ij = PETSC_TRUE;
511225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
512225daaf8SStefano Zampini       Mat T;
513225daaf8SStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,&T);CHKERRQ(ierr);
514225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
515225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
516225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
517225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
51863c07aadSStefano Zampini     }
519225daaf8SStefano Zampini   }
520225daaf8SStefano Zampini 
521225daaf8SStefano Zampini   /* we have to use hypre_Tfree to free the arrays */
52263c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
523225daaf8SStefano Zampini     void *ptrs[6] = {dii,djj,da,oii,ojj,oa};
524225daaf8SStefano Zampini     const char *names[6] = {"_hypre_csr_dii",
525225daaf8SStefano Zampini                             "_hypre_csr_djj",
526225daaf8SStefano Zampini                             "_hypre_csr_da",
527225daaf8SStefano Zampini                             "_hypre_csr_oii",
528225daaf8SStefano Zampini                             "_hypre_csr_ojj",
529225daaf8SStefano Zampini                             "_hypre_csr_oa"};
530225daaf8SStefano Zampini     for (i=0; i<6; i++) {
531225daaf8SStefano Zampini       PetscContainer c;
532225daaf8SStefano Zampini 
533225daaf8SStefano Zampini       ierr = PetscContainerCreate(comm,&c);CHKERRQ(ierr);
534225daaf8SStefano Zampini       ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr);
535225daaf8SStefano Zampini       ierr = PetscContainerSetUserDestroy(c,hypre_array_destroy);CHKERRQ(ierr);
536225daaf8SStefano Zampini       ierr = PetscObjectCompose((PetscObject)(*B),names[i],(PetscObject)c);CHKERRQ(ierr);
537225daaf8SStefano Zampini       ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
538225daaf8SStefano Zampini     }
53963c07aadSStefano Zampini   }
54063c07aadSStefano Zampini   PetscFunctionReturn(0);
54163c07aadSStefano Zampini }
54263c07aadSStefano Zampini 
54363c07aadSStefano Zampini #undef __FUNCT__
544c1a070e6SStefano Zampini #define __FUNCT__ "MatPtAP_AIJ_HYPRE"
545c1a070e6SStefano Zampini static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
546c1a070e6SStefano Zampini {
547c1a070e6SStefano Zampini   Mat_HYPRE          *hP = (Mat_HYPRE*)P->data;
548c1a070e6SStefano Zampini   hypre_ParCSRMatrix *tA,*Pparcsr,*ptapparcsr;
549c1a070e6SStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
550c1a070e6SStefano Zampini   Mat_SeqAIJ         *diag,*offd;
551c1a070e6SStefano Zampini   PetscInt           *garray,noffd,dnnz,onnz,*row_starts,*col_starts;
552c1a070e6SStefano Zampini   HYPRE_Int          type,P_owns_col_starts;
553c1a070e6SStefano Zampini   PetscBool          ismpiaij,isseqaij;
554c1a070e6SStefano Zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
555225daaf8SStefano Zampini   char               mtype[256];
556c1a070e6SStefano Zampini   PetscErrorCode     ierr;
557c1a070e6SStefano Zampini 
558c1a070e6SStefano Zampini   PetscFunctionBegin;
559c1a070e6SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
560c1a070e6SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
561c1a070e6SStefano Zampini   if (scall == MAT_REUSE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_REUSE_MATRIX");
562c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
563c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
564c1a070e6SStefano Zampini   if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type);
565225daaf8SStefano Zampini 
566225daaf8SStefano Zampini   /* It looks like we don't need to have the diagonal entries
567225daaf8SStefano Zampini      ordered first in the rows of the diagonal part
568225daaf8SStefano Zampini      for boomerAMGBuildCoarseOperator to work */
569c1a070e6SStefano Zampini   if (ismpiaij) {
570c1a070e6SStefano Zampini     Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data);
571c1a070e6SStefano Zampini 
572c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)a->A->data;
573c1a070e6SStefano Zampini     offd   = (Mat_SeqAIJ*)a->B->data;
574c1a070e6SStefano Zampini     garray = a->garray;
575c1a070e6SStefano Zampini     noffd  = a->B->cmap->N;
576c1a070e6SStefano Zampini     dnnz   = diag->nz;
577c1a070e6SStefano Zampini     onnz   = offd->nz;
578c1a070e6SStefano Zampini   } else {
579c1a070e6SStefano Zampini     diag    = (Mat_SeqAIJ*)A->data;
580c1a070e6SStefano Zampini     offd    = NULL;
581c1a070e6SStefano Zampini     garray  = NULL;
582c1a070e6SStefano Zampini     noffd   = 0;
583c1a070e6SStefano Zampini     dnnz    = diag->nz;
584c1a070e6SStefano Zampini     onnz    = 0;
585c1a070e6SStefano Zampini   }
586225daaf8SStefano Zampini 
587c1a070e6SStefano Zampini   /* create a temporary ParCSR */
588c1a070e6SStefano Zampini   if (HYPRE_AssumedPartitionCheck()) {
589c1a070e6SStefano Zampini    PetscMPIInt myid;
590c1a070e6SStefano Zampini 
591c1a070e6SStefano Zampini    ierr       = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
592c1a070e6SStefano Zampini    row_starts = A->rmap->range + myid;
593c1a070e6SStefano Zampini    col_starts = A->cmap->range + myid;
594c1a070e6SStefano Zampini   } else {
595c1a070e6SStefano Zampini    row_starts = A->rmap->range;
596c1a070e6SStefano Zampini    col_starts = A->cmap->range;
597c1a070e6SStefano Zampini   }
598c1a070e6SStefano Zampini   tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,row_starts,col_starts,noffd,dnnz,onnz);
599c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(tA,0);
600c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(tA,0);
601c1a070e6SStefano Zampini 
602225daaf8SStefano Zampini   /* set diagonal part */
603c1a070e6SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(tA);
604c1a070e6SStefano Zampini   hypre_CSRMatrixI(hdiag)           = diag->i;
605c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hdiag)           = diag->j;
606c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = diag->a;
607c1a070e6SStefano Zampini   hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz;
608c1a070e6SStefano Zampini   hypre_CSRMatrixSetRownnz(hdiag);
609c1a070e6SStefano Zampini   hypre_CSRMatrixSetDataOwner(hdiag,0);
610c1a070e6SStefano Zampini 
611225daaf8SStefano Zampini   /* set offdiagonal part */
612c1a070e6SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(tA);
613c1a070e6SStefano Zampini   if (offd) {
614c1a070e6SStefano Zampini     hypre_CSRMatrixI(hoffd)           = offd->i;
615c1a070e6SStefano Zampini     hypre_CSRMatrixJ(hoffd)           = offd->j;
616c1a070e6SStefano Zampini     hypre_CSRMatrixData(hoffd)        = offd->a;
617c1a070e6SStefano Zampini     hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz;
618c1a070e6SStefano Zampini     hypre_CSRMatrixSetRownnz(hoffd);
619c1a070e6SStefano Zampini     hypre_CSRMatrixSetDataOwner(hoffd,0);
620c1a070e6SStefano Zampini     hypre_ParCSRMatrixSetNumNonzeros(tA);
621c1a070e6SStefano Zampini     hypre_ParCSRMatrixColMapOffd(tA) = garray;
622c1a070e6SStefano Zampini   }
623c1a070e6SStefano Zampini 
624c1a070e6SStefano Zampini   /* call RAP from BoomerAMG */
625c1a070e6SStefano Zampini   /* Warning: hypre_BoomerAMGBuildCoarseOperator steals the col_starts
626c1a070e6SStefano Zampini      from Pparcsr (even if it does not own them)! */
627c1a070e6SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
628c1a070e6SStefano Zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(Pparcsr);
629c1a070e6SStefano Zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
630c1a070e6SStefano Zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(Pparcsr,tA,Pparcsr,&ptapparcsr));
631c1a070e6SStefano Zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
632c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(ptapparcsr,0);
633c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(ptapparcsr,0);
634c1a070e6SStefano Zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(Pparcsr, 1);
635c1a070e6SStefano Zampini 
636c1a070e6SStefano Zampini   /* set pointers to NULL before destroying tA */
637c1a070e6SStefano Zampini   hypre_CSRMatrixI(hdiag)          = NULL;
638c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hdiag)          = NULL;
639c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)       = NULL;
640c1a070e6SStefano Zampini   hypre_CSRMatrixI(hoffd)          = NULL;
641c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hoffd)          = NULL;
642c1a070e6SStefano Zampini   hypre_CSRMatrixData(hoffd)       = NULL;
643c1a070e6SStefano Zampini   hypre_ParCSRMatrixColMapOffd(tA) = NULL;
644c1a070e6SStefano Zampini   hypre_ParCSRMatrixDestroy(tA);
645225daaf8SStefano Zampini 
646225daaf8SStefano Zampini   /* create C depending on mtype */
647225daaf8SStefano Zampini   sprintf(mtype,MATAIJ);
648225daaf8SStefano Zampini   ierr = PetscOptionsGetString(((PetscObject)A)->options,((PetscObject)A)->prefix,"-matptap_hypre_outtype",mtype,256,NULL);CHKERRQ(ierr);
649225daaf8SStefano Zampini   ierr = MatCreateFromParCSR(ptapparcsr,mtype,PETSC_OWN_POINTER,C);CHKERRQ(ierr);
650c1a070e6SStefano Zampini   PetscFunctionReturn(0);
651c1a070e6SStefano Zampini }
652c1a070e6SStefano Zampini 
653c1a070e6SStefano Zampini #undef __FUNCT__
654cd8bc7baSStefano Zampini #define __FUNCT__ "MatPtAP_HYPRE_HYPRE"
655cd8bc7baSStefano Zampini static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
656cd8bc7baSStefano Zampini {
657cd8bc7baSStefano Zampini   hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr;
658cd8bc7baSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data, *hP = (Mat_HYPRE*)P->data;
659c1a070e6SStefano Zampini   HYPRE_Int          type,P_owns_col_starts;
660cd8bc7baSStefano Zampini   PetscErrorCode     ierr;
661cd8bc7baSStefano Zampini 
662cd8bc7baSStefano Zampini   PetscFunctionBegin;
663cd8bc7baSStefano Zampini   if (scall == MAT_REUSE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Unsupported MAT_REUSE_MATRIX");
664cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
665cd8bc7baSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
666cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
667cd8bc7baSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
668cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
669cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
670c1a070e6SStefano Zampini 
671c1a070e6SStefano Zampini   /* call RAP from BoomerAMG */
672c1a070e6SStefano Zampini   /* Warning: hypre_BoomerAMGBuildCoarseOperator steals the col_starts
673c1a070e6SStefano Zampini      from Pparcsr (even if it does not own them)! */
674c1a070e6SStefano Zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(Pparcsr);
675c1a070e6SStefano Zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
676cd8bc7baSStefano Zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr));
677cd8bc7baSStefano Zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
678c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(ptapparcsr,0);
679c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(ptapparcsr,0);
680c1a070e6SStefano Zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(Pparcsr, 1);
681c1a070e6SStefano Zampini 
682c1a070e6SStefano Zampini   /* create MatHYPRE */
683225daaf8SStefano Zampini   ierr = MatCreateFromParCSR(ptapparcsr,MATHYPRE,PETSC_OWN_POINTER,C);CHKERRQ(ierr);
684cd8bc7baSStefano Zampini   PetscFunctionReturn(0);
685cd8bc7baSStefano Zampini }
686cd8bc7baSStefano Zampini 
687cd8bc7baSStefano Zampini #undef __FUNCT__
68863c07aadSStefano Zampini #define __FUNCT__ "MatMultTranspose_HYPRE"
689ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
69063c07aadSStefano Zampini {
69163c07aadSStefano Zampini   PetscErrorCode ierr;
69263c07aadSStefano Zampini 
69363c07aadSStefano Zampini   PetscFunctionBegin;
69463c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
69563c07aadSStefano Zampini   PetscFunctionReturn(0);
69663c07aadSStefano Zampini }
69763c07aadSStefano Zampini 
69863c07aadSStefano Zampini #undef __FUNCT__
69963c07aadSStefano Zampini #define __FUNCT__ "MatMult_HYPRE"
700ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
70163c07aadSStefano Zampini {
70263c07aadSStefano Zampini   PetscErrorCode ierr;
70363c07aadSStefano Zampini 
70463c07aadSStefano Zampini   PetscFunctionBegin;
70563c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
70663c07aadSStefano Zampini   PetscFunctionReturn(0);
70763c07aadSStefano Zampini }
70863c07aadSStefano Zampini 
70963c07aadSStefano Zampini #undef __FUNCT__
71063c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_MultKernel_Private"
711ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
71263c07aadSStefano Zampini {
71363c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
71463c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
71563c07aadSStefano Zampini   hypre_ParVector    *hx,*hy;
71663c07aadSStefano Zampini   PetscScalar        *ax,*ay,*sax,*say;
71763c07aadSStefano Zampini   PetscErrorCode     ierr;
71863c07aadSStefano Zampini 
71963c07aadSStefano Zampini   PetscFunctionBegin;
72063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
72163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
72263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
72363c07aadSStefano Zampini   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
72463c07aadSStefano Zampini   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
72563c07aadSStefano Zampini   if (trans) {
72658968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ay,say);
72758968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax);
72863c07aadSStefano Zampini     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
72958968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,say,ay);
73058968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax);
73163c07aadSStefano Zampini   } else {
73258968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax);
73358968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ay,say);
73463c07aadSStefano Zampini     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
73558968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax);
73658968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,say,ay);
73763c07aadSStefano Zampini   }
73863c07aadSStefano Zampini   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
73963c07aadSStefano Zampini   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
74063c07aadSStefano Zampini   PetscFunctionReturn(0);
74163c07aadSStefano Zampini }
74263c07aadSStefano Zampini 
74363c07aadSStefano Zampini #undef __FUNCT__
74463c07aadSStefano Zampini #define __FUNCT__ "MatDestroy_HYPRE"
745ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A)
74663c07aadSStefano Zampini {
74763c07aadSStefano Zampini   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
74863c07aadSStefano Zampini   PetscErrorCode ierr;
74963c07aadSStefano Zampini 
75063c07aadSStefano Zampini   PetscFunctionBegin;
75163c07aadSStefano Zampini   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
75263c07aadSStefano Zampini   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
753978814f1SStefano Zampini   if (hA->ij) {
754978814f1SStefano Zampini     if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL;
755978814f1SStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
756978814f1SStefano Zampini   }
75763c07aadSStefano Zampini   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); }
75863c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
7592df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr);
760c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr);
761c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr);
762d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPRESetPreallocation_C",NULL);CHKERRQ(ierr);
763*dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPREGetParCSR_C",NULL);CHKERRQ(ierr);
76463c07aadSStefano Zampini   ierr = PetscFree(A->data);CHKERRQ(ierr);
76563c07aadSStefano Zampini   PetscFunctionReturn(0);
76663c07aadSStefano Zampini }
76763c07aadSStefano Zampini 
76863c07aadSStefano Zampini #undef __FUNCT__
76963c07aadSStefano Zampini #define __FUNCT__ "MatSetUp_HYPRE"
770ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A)
77163c07aadSStefano Zampini {
77263c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
77363c07aadSStefano Zampini   Vec                x,b;
77463c07aadSStefano Zampini   PetscErrorCode     ierr;
77563c07aadSStefano Zampini 
77663c07aadSStefano Zampini   PetscFunctionBegin;
77763c07aadSStefano Zampini   if (hA->x) PetscFunctionReturn(0);
77863c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
77963c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
78063c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
78163c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
78263c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
78363c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
78463c07aadSStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
78563c07aadSStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
78663c07aadSStefano Zampini   PetscFunctionReturn(0);
78763c07aadSStefano Zampini }
78863c07aadSStefano Zampini 
78963c07aadSStefano Zampini #undef __FUNCT__
79063c07aadSStefano Zampini #define __FUNCT__ "MatAssemblyEnd_HYPRE"
791ea9daf28SStefano Zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
79263c07aadSStefano Zampini {
79363c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
79463c07aadSStefano Zampini   PetscErrorCode     ierr;
79563c07aadSStefano Zampini 
79663c07aadSStefano Zampini   PetscFunctionBegin;
797d975228cSstefano_zampini   if (mode == MAT_FLUSH_ASSEMBLY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_FLUSH_ASSEMBLY currently not supported with MATHYPRE");
79863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
79963c07aadSStefano Zampini   PetscFunctionReturn(0);
80063c07aadSStefano Zampini }
80163c07aadSStefano Zampini 
802d975228cSstefano_zampini #define MATHYPRE_SCRATCH 2048
803d975228cSstefano_zampini 
804d975228cSstefano_zampini #undef __FUNCT__
805d975228cSstefano_zampini #define __FUNCT__ "MatSetValues_HYPRE"
806d975228cSstefano_zampini PetscErrorCode MatSetValues_HYPRE(Mat A, PetscInt nr, const PetscInt rows[], PetscInt nc, const PetscInt cols[], const PetscScalar v[], InsertMode ins)
807d975228cSstefano_zampini {
808d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
809d975228cSstefano_zampini   PetscScalar        *vals = (PetscScalar *)v;
810d975228cSstefano_zampini   PetscScalar        sscr[MATHYPRE_SCRATCH];
811d975228cSstefano_zampini   PetscInt           cscr[2][MATHYPRE_SCRATCH];
812d975228cSstefano_zampini   PetscInt           i,nzc;
813d975228cSstefano_zampini   PetscErrorCode     ierr;
814d975228cSstefano_zampini 
815d975228cSstefano_zampini   PetscFunctionBegin;
816d975228cSstefano_zampini   for (i=0,nzc=0;i<nc;i++) {
817d975228cSstefano_zampini     if (cols[i] >= 0) {
818d975228cSstefano_zampini       cscr[0][nzc  ] = cols[i];
819d975228cSstefano_zampini       cscr[1][nzc++] = i;
820d975228cSstefano_zampini     }
821d975228cSstefano_zampini   }
822d975228cSstefano_zampini   if (!nzc) PetscFunctionReturn(0);
823d975228cSstefano_zampini 
824d975228cSstefano_zampini   if (ins == ADD_VALUES) {
825d975228cSstefano_zampini     for (i=0;i<nr;i++) {
826d975228cSstefano_zampini       if (rows[i] >= 0) {
827d975228cSstefano_zampini         PetscInt j;
828d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
829d975228cSstefano_zampini         PetscStackCallStandard(HYPRE_IJMatrixAddToValues,(hA->ij,1,&nzc,rows+i,cscr[0],sscr));
830d975228cSstefano_zampini       }
831d975228cSstefano_zampini       vals += nc;
832d975228cSstefano_zampini     }
833d975228cSstefano_zampini   } else { /* INSERT_VALUES */
834d975228cSstefano_zampini #if defined(PETSC_USE_DEBUG)
835d975228cSstefano_zampini     /* Insert values cannot be used to insert offproc entries */
836d975228cSstefano_zampini     PetscInt rst,ren;
837d975228cSstefano_zampini     ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
838d975228cSstefano_zampini     for (i=0;i<nr;i++)
839d975228cSstefano_zampini       if (rows[i] >= 0 && (rows[i] < rst || rows[i] >= ren)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use INSERT_VALUES for off-proc entries with MatHYPRE. Use ADD_VALUES instead");
840d975228cSstefano_zampini #endif
841d975228cSstefano_zampini     for (i=0;i<nr;i++) {
842d975228cSstefano_zampini       if (rows[i] >= 0) {
843d975228cSstefano_zampini         PetscInt j;
844d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
845d975228cSstefano_zampini         PetscStackCallStandard(HYPRE_IJMatrixSetValues,(hA->ij,1,&nzc,rows+i,cscr[0],sscr));
846d975228cSstefano_zampini       }
847d975228cSstefano_zampini       vals += nc;
848d975228cSstefano_zampini     }
849d975228cSstefano_zampini   }
850d975228cSstefano_zampini   PetscFunctionReturn(0);
851d975228cSstefano_zampini }
852d975228cSstefano_zampini 
853d975228cSstefano_zampini #undef __FUNCT__
854d975228cSstefano_zampini #define __FUNCT__ "MatHYPRESetPreallocation_HYPRE"
855d975228cSstefano_zampini static PetscErrorCode MatHYPRESetPreallocation_HYPRE(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
856d975228cSstefano_zampini {
857d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
858d975228cSstefano_zampini   PetscInt           i,*hdnnz,*honnz;
859d975228cSstefano_zampini   PetscInt           rs,re,cs,ce;
860d975228cSstefano_zampini   PetscMPIInt        size;
861d975228cSstefano_zampini   PetscErrorCode     ierr;
862d975228cSstefano_zampini 
863d975228cSstefano_zampini   PetscFunctionBegin;
864d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
865d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
866d975228cSstefano_zampini   rs   = A->rmap->rstart;
867d975228cSstefano_zampini   re   = A->rmap->rend;
868d975228cSstefano_zampini   cs   = A->cmap->rstart;
869d975228cSstefano_zampini   ce   = A->cmap->rend;
870d975228cSstefano_zampini   if (!hA->ij) {
871d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rs,re-1,cs,ce-1,&hA->ij));
872d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
873d975228cSstefano_zampini   } else {
874d975228cSstefano_zampini     PetscInt hrs,hre,hcs,hce;
875d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixGetLocalRange,(hA->ij,&hrs,&hre,&hcs,&hce));
876d975228cSstefano_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);
877d975228cSstefano_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);
878d975228cSstefano_zampini   }
879d975228cSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
880d975228cSstefano_zampini 
881d975228cSstefano_zampini   if (!dnnz) {
882d975228cSstefano_zampini     ierr = PetscMalloc1(A->rmap->n,&hdnnz);CHKERRQ(ierr);
883d975228cSstefano_zampini     for (i=0;i<A->rmap->n;i++) hdnnz[i] = dnz;
884d975228cSstefano_zampini   } else {
885d975228cSstefano_zampini     hdnnz = (PetscInt*)dnnz;
886d975228cSstefano_zampini   }
887d975228cSstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
888d975228cSstefano_zampini   if (size > 1) {
889d975228cSstefano_zampini     if (!onnz) {
890d975228cSstefano_zampini       ierr = PetscMalloc1(A->rmap->n,&honnz);CHKERRQ(ierr);
891d975228cSstefano_zampini       for (i=0;i<A->rmap->n;i++) honnz[i] = onz;
892d975228cSstefano_zampini     } else {
893d975228cSstefano_zampini       honnz = (PetscInt*)onnz;
894d975228cSstefano_zampini     }
895d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(hA->ij,hdnnz,honnz));
896d975228cSstefano_zampini   } else {
897d975228cSstefano_zampini     honnz = NULL;
898d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetRowSizes,(hA->ij,hdnnz));
899d975228cSstefano_zampini   }
900d975228cSstefano_zampini   if (!dnnz) {
901d975228cSstefano_zampini     ierr = PetscFree(hdnnz);CHKERRQ(ierr);
902d975228cSstefano_zampini   }
903d975228cSstefano_zampini   if (!onnz && honnz) {
904d975228cSstefano_zampini     ierr = PetscFree(honnz);CHKERRQ(ierr);
905d975228cSstefano_zampini   }
906d975228cSstefano_zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
907d975228cSstefano_zampini 
908d975228cSstefano_zampini   /* SetDiagOffdSizes sets hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0 */
909d975228cSstefano_zampini   {
910d975228cSstefano_zampini     hypre_AuxParCSRMatrix *aux_matrix;
911d975228cSstefano_zampini     aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij);
912d975228cSstefano_zampini     hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1;
913d975228cSstefano_zampini   }
914d975228cSstefano_zampini   PetscFunctionReturn(0);
915d975228cSstefano_zampini }
916d975228cSstefano_zampini 
917d975228cSstefano_zampini /*@C
918d975228cSstefano_zampini    MatHYPRESetPreallocation - Preallocates memory for a sparse parallel matrix in HYPRE IJ format
919d975228cSstefano_zampini 
920d975228cSstefano_zampini    Collective on Mat
921d975228cSstefano_zampini 
922d975228cSstefano_zampini    Input Parameters:
923d975228cSstefano_zampini +  A - the matrix
924d975228cSstefano_zampini .  dnz  - number of nonzeros per row in DIAGONAL portion of local submatrix
925d975228cSstefano_zampini           (same value is used for all local rows)
926d975228cSstefano_zampini .  dnnz - array containing the number of nonzeros in the various rows of the
927d975228cSstefano_zampini           DIAGONAL portion of the local submatrix (possibly different for each row)
928d975228cSstefano_zampini           or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure.
929d975228cSstefano_zampini           The size of this array is equal to the number of local rows, i.e 'm'.
930d975228cSstefano_zampini           For matrices that will be factored, you must leave room for (and set)
931d975228cSstefano_zampini           the diagonal entry even if it is zero.
932d975228cSstefano_zampini .  onz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
933d975228cSstefano_zampini           submatrix (same value is used for all local rows).
934d975228cSstefano_zampini -  onnz - array containing the number of nonzeros in the various rows of the
935d975228cSstefano_zampini           OFF-DIAGONAL portion of the local submatrix (possibly different for
936d975228cSstefano_zampini           each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero
937d975228cSstefano_zampini           structure. The size of this array is equal to the number
938d975228cSstefano_zampini           of local rows, i.e 'm'.
939d975228cSstefano_zampini 
940d975228cSstefano_zampini    Notes: If the *nnz parameter is given then the *nz parameter is ignored; for sequential matrices, onz and onnz are ignored.
941d975228cSstefano_zampini 
942d975228cSstefano_zampini    Level: intermediate
943d975228cSstefano_zampini 
944d975228cSstefano_zampini .keywords: matrix, aij, compressed row, sparse, parallel
945d975228cSstefano_zampini 
946d975228cSstefano_zampini .seealso: MatCreate(), MatMPIAIJSetPreallocation, MATHYPRE
947d975228cSstefano_zampini @*/
948d975228cSstefano_zampini #undef __FUNCT__
949d975228cSstefano_zampini #define __FUNCT__ "MatHYPRESetPreallocation"
950d975228cSstefano_zampini PetscErrorCode MatHYPRESetPreallocation(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
951d975228cSstefano_zampini {
952d975228cSstefano_zampini   PetscErrorCode ierr;
953d975228cSstefano_zampini 
954d975228cSstefano_zampini   PetscFunctionBegin;
955d975228cSstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
956d975228cSstefano_zampini   PetscValidType(A,1);
957d975228cSstefano_zampini   ierr = PetscTryMethod(A,"MatHYPRESetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(A,dnz,dnnz,onz,onnz));CHKERRQ(ierr);
958d975228cSstefano_zampini   PetscFunctionReturn(0);
959d975228cSstefano_zampini }
960d975228cSstefano_zampini 
961225daaf8SStefano Zampini /*
962225daaf8SStefano Zampini    MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix
963225daaf8SStefano Zampini 
964225daaf8SStefano Zampini    Collective
965225daaf8SStefano Zampini 
966225daaf8SStefano Zampini    Input Parameters:
967225daaf8SStefano Zampini +  vparcsr  - the pointer to the hypre_ParCSRMatrix
968bb4689ddSStefano Zampini .  mtype    - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported.
969225daaf8SStefano Zampini -  copymode - PETSc copying options
970225daaf8SStefano Zampini 
971225daaf8SStefano Zampini    Output Parameter:
972225daaf8SStefano Zampini .  A  - the matrix
973225daaf8SStefano Zampini 
974225daaf8SStefano Zampini    Level: intermediate
975225daaf8SStefano Zampini 
976225daaf8SStefano Zampini .seealso: MatHYPRE, PetscCopyMode
977225daaf8SStefano Zampini */
978978814f1SStefano Zampini #undef __FUNCT__
979225daaf8SStefano Zampini #define __FUNCT__ "MatCreateFromParCSR"
980*dd9c0a25Sstefano_zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(hypre_ParCSRMatrix *vparcsr, MatType mtype, PetscCopyMode copymode, Mat* A)
981978814f1SStefano Zampini {
982225daaf8SStefano Zampini   Mat                   T;
983978814f1SStefano Zampini   Mat_HYPRE             *hA;
98463c07aadSStefano Zampini   hypre_ParCSRMatrix    *parcsr;
985978814f1SStefano Zampini   MPI_Comm              comm;
986978814f1SStefano Zampini   PetscInt              rstart,rend,cstart,cend,M,N;
987bb4689ddSStefano Zampini   PetscBool             isseqaij,ismpiaij,isaij,ishyp,isis;
988978814f1SStefano Zampini   PetscErrorCode        ierr;
989978814f1SStefano Zampini 
990978814f1SStefano Zampini   PetscFunctionBegin;
991978814f1SStefano Zampini   parcsr = (hypre_ParCSRMatrix *)vparcsr;
992978814f1SStefano Zampini   comm   = hypre_ParCSRMatrixComm(parcsr);
993225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
994225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
995225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr);
996225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr);
9978cfe8d00SStefano Zampini   ierr   = PetscStrcmp(mtype,MATIS,&isis);CHKERRQ(ierr);
998225daaf8SStefano Zampini   isaij  = (PetscBool)(isseqaij || ismpiaij || isaij);
999bb4689ddSStefano 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);
1000225daaf8SStefano Zampini   if (ishyp && copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES");
1001978814f1SStefano Zampini 
1002978814f1SStefano Zampini   /* access ParCSRMatrix */
1003978814f1SStefano Zampini   rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr);
1004978814f1SStefano Zampini   rend   = hypre_ParCSRMatrixLastRowIndex(parcsr);
1005978814f1SStefano Zampini   cstart = hypre_ParCSRMatrixFirstColDiag(parcsr);
1006978814f1SStefano Zampini   cend   = hypre_ParCSRMatrixLastColDiag(parcsr);
1007978814f1SStefano Zampini   M      = hypre_ParCSRMatrixGlobalNumRows(parcsr);
1008978814f1SStefano Zampini   N      = hypre_ParCSRMatrixGlobalNumCols(parcsr);
1009978814f1SStefano Zampini 
1010978814f1SStefano Zampini   /* create PETSc matrix with MatHYPRE */
1011225daaf8SStefano Zampini   ierr = MatCreate(comm,&T);CHKERRQ(ierr);
1012225daaf8SStefano Zampini   ierr = MatSetSizes(T,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr);
1013225daaf8SStefano Zampini   ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr);
1014225daaf8SStefano Zampini   ierr = MatSetUp(T);CHKERRQ(ierr);
1015225daaf8SStefano Zampini   hA   = (Mat_HYPRE*)(T->data);
1016978814f1SStefano Zampini 
1017978814f1SStefano Zampini   /* create HYPRE_IJMatrix */
1018978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
1019978814f1SStefano Zampini 
1020978814f1SStefano Zampini   /* set ParCSR object */
1021978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1022978814f1SStefano Zampini   hypre_IJMatrixObject(hA->ij) = parcsr;
1023978814f1SStefano Zampini 
1024978814f1SStefano Zampini   /* set assembled flag */
1025978814f1SStefano Zampini   hypre_IJMatrixAssembleFlag(hA->ij) = 1;
1026978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1027225daaf8SStefano Zampini   if (ishyp) {
1028978814f1SStefano Zampini     /* prevent from freeing the pointer */
1029978814f1SStefano Zampini     if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE;
1030225daaf8SStefano Zampini     *A = T;
1031978814f1SStefano Zampini     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1032978814f1SStefano Zampini     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1033bb4689ddSStefano Zampini   } else if (isaij) {
1034bb4689ddSStefano Zampini     if (copymode != PETSC_OWN_POINTER) {
1035225daaf8SStefano Zampini       /* prevent from freeing the pointer */
1036225daaf8SStefano Zampini       hA->inner_free = PETSC_FALSE;
1037225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
1038225daaf8SStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
1039225daaf8SStefano Zampini     } else { /* AIJ return type with PETSC_OWN_POINTER */
1040225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1041225daaf8SStefano Zampini       *A   = T;
1042225daaf8SStefano Zampini     }
1043bb4689ddSStefano Zampini   } else if (isis) {
10448cfe8d00SStefano Zampini     ierr = MatConvert_HYPRE_IS(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
10458cfe8d00SStefano Zampini     if (copymode != PETSC_OWN_POINTER) hA->inner_free = PETSC_FALSE;
1046bb4689ddSStefano Zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
1047bb4689ddSStefano Zampini   }
1048978814f1SStefano Zampini   PetscFunctionReturn(0);
1049978814f1SStefano Zampini }
1050978814f1SStefano Zampini 
105163c07aadSStefano Zampini #undef __FUNCT__
1052*dd9c0a25Sstefano_zampini #define __FUNCT__ "MatHYPREGetParCSR_HYPRE"
1053*dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR_HYPRE(Mat A, hypre_ParCSRMatrix **parcsr)
1054*dd9c0a25Sstefano_zampini {
1055*dd9c0a25Sstefano_zampini   Mat_HYPRE*            hA = (Mat_HYPRE*)A->data;
1056*dd9c0a25Sstefano_zampini   HYPRE_Int             type;
1057*dd9c0a25Sstefano_zampini   PetscErrorCode        ierr;
1058*dd9c0a25Sstefano_zampini 
1059*dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1060*dd9c0a25Sstefano_zampini   if (!hA->ij) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"HYPRE_IJMatrix not present");
1061*dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
1062*dd9c0a25Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"HYPRE_IJMatrix is not of type HYPRE_PARCSR");
1063*dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)parcsr));
1064*dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1065*dd9c0a25Sstefano_zampini }
1066*dd9c0a25Sstefano_zampini 
1067*dd9c0a25Sstefano_zampini /*
1068*dd9c0a25Sstefano_zampini    MatHYPREGetParCSR - Gets the pointer to the ParCSR matrix
1069*dd9c0a25Sstefano_zampini 
1070*dd9c0a25Sstefano_zampini    Not collective
1071*dd9c0a25Sstefano_zampini 
1072*dd9c0a25Sstefano_zampini    Input Parameters:
1073*dd9c0a25Sstefano_zampini +  A  - the MATHYPRE object
1074*dd9c0a25Sstefano_zampini 
1075*dd9c0a25Sstefano_zampini    Output Parameter:
1076*dd9c0a25Sstefano_zampini .  parcsr  - the pointer to the hypre_ParCSRMatrix
1077*dd9c0a25Sstefano_zampini 
1078*dd9c0a25Sstefano_zampini    Level: intermediate
1079*dd9c0a25Sstefano_zampini 
1080*dd9c0a25Sstefano_zampini .seealso: MatHYPRE, PetscCopyMode
1081*dd9c0a25Sstefano_zampini */
1082*dd9c0a25Sstefano_zampini #undef __FUNCT__
1083*dd9c0a25Sstefano_zampini #define __FUNCT__ "MatHYPREGetParCSR"
1084*dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR(Mat A, hypre_ParCSRMatrix **parcsr)
1085*dd9c0a25Sstefano_zampini {
1086*dd9c0a25Sstefano_zampini   PetscErrorCode ierr;
1087*dd9c0a25Sstefano_zampini 
1088*dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1089*dd9c0a25Sstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1090*dd9c0a25Sstefano_zampini   PetscValidType(A,1);
1091*dd9c0a25Sstefano_zampini   ierr = PetscUseMethod(A,"MatHYPREGetParCSR_C",(Mat,hypre_ParCSRMatrix**),(A,parcsr));CHKERRQ(ierr);
1092*dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1093*dd9c0a25Sstefano_zampini }
1094*dd9c0a25Sstefano_zampini 
1095*dd9c0a25Sstefano_zampini #undef __FUNCT__
109663c07aadSStefano Zampini #define __FUNCT__ "MatCreate_HYPRE"
109763c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
109863c07aadSStefano Zampini {
109963c07aadSStefano Zampini   Mat_HYPRE      *hB;
110063c07aadSStefano Zampini   PetscErrorCode ierr;
110163c07aadSStefano Zampini 
110263c07aadSStefano Zampini   PetscFunctionBegin;
110363c07aadSStefano Zampini   ierr           = PetscNewLog(B,&hB);CHKERRQ(ierr);
1104978814f1SStefano Zampini   hB->inner_free = PETSC_TRUE;
1105978814f1SStefano Zampini 
110663c07aadSStefano Zampini   B->data       = (void*)hB;
110763c07aadSStefano Zampini   B->rmap->bs   = 1;
110863c07aadSStefano Zampini   B->assembled  = PETSC_FALSE;
110963c07aadSStefano Zampini 
111063c07aadSStefano Zampini   B->ops->mult          = MatMult_HYPRE;
111163c07aadSStefano Zampini   B->ops->multtranspose = MatMultTranspose_HYPRE;
111263c07aadSStefano Zampini   B->ops->setup         = MatSetUp_HYPRE;
111363c07aadSStefano Zampini   B->ops->destroy       = MatDestroy_HYPRE;
111463c07aadSStefano Zampini   B->ops->assemblyend   = MatAssemblyEnd_HYPRE;
1115cd8bc7baSStefano Zampini   B->ops->ptap          = MatPtAP_HYPRE_HYPRE;
1116d975228cSstefano_zampini   B->ops->setvalues     = MatSetValues_HYPRE;
111763c07aadSStefano Zampini 
111863c07aadSStefano Zampini   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
111963c07aadSStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
112063c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
11212df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr);
1122c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1123c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1124d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPRESetPreallocation_C",MatHYPRESetPreallocation_HYPRE);CHKERRQ(ierr);
1125*dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPREGetParCSR_C",MatHYPREGetParCSR_HYPRE);CHKERRQ(ierr);
112663c07aadSStefano Zampini   PetscFunctionReturn(0);
112763c07aadSStefano Zampini }
112863c07aadSStefano Zampini 
1129225daaf8SStefano Zampini #undef __FUNCT__
1130225daaf8SStefano Zampini #define __FUNCT__ "hypre_array_destroy"
1131225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void *ptr)
1132225daaf8SStefano Zampini {
1133225daaf8SStefano Zampini    PetscFunctionBegin;
1134225daaf8SStefano Zampini    hypre_TFree(ptr);
1135225daaf8SStefano Zampini    PetscFunctionReturn(0);
1136225daaf8SStefano Zampini }
1137225daaf8SStefano Zampini 
113863c07aadSStefano Zampini #if 0
113963c07aadSStefano Zampini /*
114063c07aadSStefano Zampini     Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format
114163c07aadSStefano Zampini 
114263c07aadSStefano Zampini     This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first
114363c07aadSStefano Zampini     which will corrupt the PETSc data structure if we did this. Need a work around to this problem.
114463c07aadSStefano Zampini */
114563c07aadSStefano Zampini #include <_hypre_IJ_mv.h>
114663c07aadSStefano Zampini #include <HYPRE_IJ_mv.h>
114763c07aadSStefano Zampini 
114863c07aadSStefano Zampini #undef __FUNCT__
114963c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixLink"
115063c07aadSStefano Zampini PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij)
115163c07aadSStefano Zampini {
115263c07aadSStefano Zampini   PetscErrorCode        ierr;
115363c07aadSStefano Zampini   PetscInt              rstart,rend,cstart,cend;
115463c07aadSStefano Zampini   PetscBool             flg;
115563c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
115663c07aadSStefano Zampini 
115763c07aadSStefano Zampini   PetscFunctionBegin;
115863c07aadSStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
115963c07aadSStefano Zampini   PetscValidType(A,1);
116063c07aadSStefano Zampini   PetscValidPointer(ij,2);
116163c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
116263c07aadSStefano Zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices");
116363c07aadSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
116463c07aadSStefano Zampini 
116563c07aadSStefano Zampini   ierr   = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
116663c07aadSStefano Zampini   rstart = A->rmap->rstart;
116763c07aadSStefano Zampini   rend   = A->rmap->rend;
116863c07aadSStefano Zampini   cstart = A->cmap->rstart;
116963c07aadSStefano Zampini   cend   = A->cmap->rend;
117063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij));
117163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
117263c07aadSStefano Zampini 
117363c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij));
117463c07aadSStefano Zampini   PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij));
117563c07aadSStefano Zampini 
117663c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
117763c07aadSStefano Zampini 
117863c07aadSStefano Zampini   /* this is the Hack part where we monkey directly with the hypre datastructures */
117963c07aadSStefano Zampini 
118063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij));
118163c07aadSStefano Zampini   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
118263c07aadSStefano Zampini   PetscFunctionReturn(0);
118363c07aadSStefano Zampini }
118463c07aadSStefano Zampini #endif
1185