xref: /petsc/src/mat/impls/hypre/mhypre.c (revision 3dad0653c8b7d79802f38d512f3aec944c073f92)
163c07aadSStefano Zampini 
263c07aadSStefano Zampini /*
363c07aadSStefano Zampini     Creates hypre ijmatrix from PETSc matrix
463c07aadSStefano Zampini */
5225daaf8SStefano Zampini 
6225daaf8SStefano Zampini /*MC
7225daaf8SStefano Zampini    MATHYPRE - MATHYPRE = "hypre" - A matrix type to be used for sequential and parallel sparse matrices
8225daaf8SStefano Zampini           based on the hypre IJ interface.
9225daaf8SStefano Zampini 
10225daaf8SStefano Zampini    Level: intermediate
11225daaf8SStefano Zampini 
12225daaf8SStefano Zampini .seealso: MatCreate()
13225daaf8SStefano Zampini M*/
14225daaf8SStefano Zampini 
15dd9c0a25Sstefano_zampini #include <petscmathypre.h>
1663c07aadSStefano Zampini #include <petsc/private/matimpl.h>
1763c07aadSStefano Zampini #include <../src/mat/impls/hypre/mhypre.h>
1863c07aadSStefano Zampini #include <../src/mat/impls/aij/mpi/mpiaij.h>
1958968eb6SStefano Zampini #include <../src/vec/vec/impls/hypre/vhyp.h>
2058968eb6SStefano Zampini #include <HYPRE.h>
21c1a070e6SStefano Zampini #include <HYPRE_utilities.h>
22cd8bc7baSStefano Zampini #include <_hypre_parcsr_ls.h>
2363c07aadSStefano Zampini 
2463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*);
2563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix);
2663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
2763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
2863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool);
29225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void*);
3063c07aadSStefano Zampini 
3163c07aadSStefano Zampini #undef __FUNCT__
3263c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate"
3363c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij)
3463c07aadSStefano Zampini {
3563c07aadSStefano Zampini   PetscErrorCode ierr;
3663c07aadSStefano Zampini   PetscInt       i,n_d,n_o;
3763c07aadSStefano Zampini   const PetscInt *ia_d,*ia_o;
3863c07aadSStefano Zampini   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
3963c07aadSStefano Zampini   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
4063c07aadSStefano Zampini 
4163c07aadSStefano Zampini   PetscFunctionBegin;
4263c07aadSStefano Zampini   if (A_d) { /* determine number of nonzero entries in local diagonal part */
4363c07aadSStefano Zampini     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
4463c07aadSStefano Zampini     if (done_d) {
4563c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr);
4663c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
4763c07aadSStefano Zampini         nnz_d[i] = ia_d[i+1] - ia_d[i];
4863c07aadSStefano Zampini       }
4963c07aadSStefano Zampini     }
5063c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
5163c07aadSStefano Zampini   }
5263c07aadSStefano Zampini   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
5363c07aadSStefano Zampini     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
5463c07aadSStefano Zampini     if (done_o) {
5563c07aadSStefano Zampini       ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr);
5663c07aadSStefano Zampini       for (i=0; i<n_o; i++) {
5763c07aadSStefano Zampini         nnz_o[i] = ia_o[i+1] - ia_o[i];
5863c07aadSStefano Zampini       }
5963c07aadSStefano Zampini     }
6063c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
6163c07aadSStefano Zampini   }
6263c07aadSStefano Zampini   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
6363c07aadSStefano Zampini     if (!done_o) { /* only diagonal part */
6463c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
6563c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
6663c07aadSStefano Zampini         nnz_o[i] = 0;
6763c07aadSStefano Zampini       }
6863c07aadSStefano Zampini     }
6963c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o));
7063c07aadSStefano Zampini     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
7163c07aadSStefano Zampini     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
7263c07aadSStefano Zampini   }
7363c07aadSStefano Zampini   PetscFunctionReturn(0);
7463c07aadSStefano Zampini }
7563c07aadSStefano Zampini 
7663c07aadSStefano Zampini #undef __FUNCT__
7763c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_CreateFromMat"
7863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA)
7963c07aadSStefano Zampini {
8063c07aadSStefano Zampini   PetscErrorCode ierr;
8163c07aadSStefano Zampini   PetscInt       rstart,rend,cstart,cend;
8263c07aadSStefano Zampini 
8363c07aadSStefano Zampini   PetscFunctionBegin;
84d975228cSstefano_zampini   ierr   = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
85d975228cSstefano_zampini   ierr   = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
8663c07aadSStefano Zampini   rstart = A->rmap->rstart;
8763c07aadSStefano Zampini   rend   = A->rmap->rend;
8863c07aadSStefano Zampini   cstart = A->cmap->rstart;
8963c07aadSStefano Zampini   cend   = A->cmap->rend;
9063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
9163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
9263c07aadSStefano Zampini   {
9363c07aadSStefano Zampini     PetscBool      same;
9463c07aadSStefano Zampini     Mat            A_d,A_o;
9563c07aadSStefano Zampini     const PetscInt *colmap;
9663c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
9763c07aadSStefano Zampini     if (same) {
9863c07aadSStefano Zampini       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
9963c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
10063c07aadSStefano Zampini       PetscFunctionReturn(0);
10163c07aadSStefano Zampini     }
10263c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
10363c07aadSStefano Zampini     if (same) {
10463c07aadSStefano Zampini       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
10563c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
10663c07aadSStefano Zampini       PetscFunctionReturn(0);
10763c07aadSStefano Zampini     }
10863c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
10963c07aadSStefano Zampini     if (same) {
11063c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
11163c07aadSStefano Zampini       PetscFunctionReturn(0);
11263c07aadSStefano Zampini     }
11363c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
11463c07aadSStefano Zampini     if (same) {
11563c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
11663c07aadSStefano Zampini       PetscFunctionReturn(0);
11763c07aadSStefano Zampini     }
11863c07aadSStefano Zampini   }
11963c07aadSStefano Zampini   PetscFunctionReturn(0);
12063c07aadSStefano Zampini }
12163c07aadSStefano Zampini 
12263c07aadSStefano Zampini #undef __FUNCT__
12363c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixCopy"
12463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij)
12563c07aadSStefano Zampini {
12663c07aadSStefano Zampini   PetscErrorCode    ierr;
12763c07aadSStefano Zampini   PetscInt          i,rstart,rend,ncols,nr,nc;
12863c07aadSStefano Zampini   const PetscScalar *values;
12963c07aadSStefano Zampini   const PetscInt    *cols;
13063c07aadSStefano Zampini   PetscBool         flg;
13163c07aadSStefano Zampini 
13263c07aadSStefano Zampini   PetscFunctionBegin;
13363c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
13463c07aadSStefano Zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
13563c07aadSStefano Zampini   if (flg && nr == nc) {
13663c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
13763c07aadSStefano Zampini     PetscFunctionReturn(0);
13863c07aadSStefano Zampini   }
13963c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
14063c07aadSStefano Zampini   if (flg) {
14163c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
14263c07aadSStefano Zampini     PetscFunctionReturn(0);
14363c07aadSStefano Zampini   }
14463c07aadSStefano Zampini 
14563c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
14663c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
14763c07aadSStefano Zampini   for (i=rstart; i<rend; i++) {
14863c07aadSStefano Zampini     ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
149e3977e59Sstefano_zampini     if (ncols) {
15063c07aadSStefano Zampini       PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
151e3977e59Sstefano_zampini     }
15263c07aadSStefano Zampini     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
15363c07aadSStefano Zampini   }
15463c07aadSStefano Zampini   PetscFunctionReturn(0);
15563c07aadSStefano Zampini }
15663c07aadSStefano Zampini 
15763c07aadSStefano Zampini #undef __FUNCT__
15863c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ"
15963c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij)
16063c07aadSStefano Zampini {
16163c07aadSStefano Zampini   PetscErrorCode        ierr;
16263c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag = (Mat_SeqAIJ*)A->data;
16358968eb6SStefano Zampini   HYPRE_Int             type;
16463c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
16563c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
16663c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag;
16763c07aadSStefano Zampini 
16863c07aadSStefano Zampini   PetscFunctionBegin;
16963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
170ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
171ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
172ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
17363c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
17463c07aadSStefano Zampini   /*
17563c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
17663c07aadSStefano Zampini   */
17763c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));
17863c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));
17963c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
180ea9daf28SStefano Zampini 
181ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
18263c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
18363c07aadSStefano Zampini   PetscFunctionReturn(0);
18463c07aadSStefano Zampini }
18563c07aadSStefano Zampini 
18663c07aadSStefano Zampini #undef __FUNCT__
18763c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ"
18863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
18963c07aadSStefano Zampini {
19063c07aadSStefano Zampini   PetscErrorCode        ierr;
19163c07aadSStefano Zampini   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
19263c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag,*poffd;
19363c07aadSStefano Zampini   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
19458968eb6SStefano Zampini   HYPRE_Int             type;
19563c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
19663c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
19763c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag,*hoffd;
19863c07aadSStefano Zampini 
19963c07aadSStefano Zampini   PetscFunctionBegin;
20063c07aadSStefano Zampini   pdiag = (Mat_SeqAIJ*) pA->A->data;
20163c07aadSStefano Zampini   poffd = (Mat_SeqAIJ*) pA->B->data;
20263c07aadSStefano Zampini   /* cstart is only valid for square MPIAIJ layed out in the usual way */
20363c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
20463c07aadSStefano Zampini 
20563c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
206ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
207ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
208ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
20963c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
21063c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
21163c07aadSStefano Zampini 
21263c07aadSStefano Zampini   /*
21363c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
21463c07aadSStefano Zampini   */
21563c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
21663c07aadSStefano Zampini   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
21763c07aadSStefano Zampini   jj  = (PetscInt*)hdiag->j;
21863c07aadSStefano Zampini   pjj = (PetscInt*)pdiag->j;
21963c07aadSStefano Zampini   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
22063c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
22163c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
22263c07aadSStefano Zampini   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
22363c07aadSStefano Zampini      If we hacked a hypre a bit more we might be able to avoid this step */
22463c07aadSStefano Zampini   jj  = (PetscInt*) hoffd->j;
22563c07aadSStefano Zampini   pjj = (PetscInt*) poffd->j;
22663c07aadSStefano Zampini   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
22763c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
22863c07aadSStefano Zampini 
229ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
23063c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
23163c07aadSStefano Zampini   PetscFunctionReturn(0);
23263c07aadSStefano Zampini }
23363c07aadSStefano Zampini 
23463c07aadSStefano Zampini #undef __FUNCT__
2352df22349SStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_IS"
2362df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B)
2372df22349SStefano Zampini {
2382df22349SStefano Zampini   Mat_HYPRE*             mhA = (Mat_HYPRE*)(A->data);
2392df22349SStefano Zampini   Mat                    lA;
2402df22349SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
2412df22349SStefano Zampini   IS                     is;
2422df22349SStefano Zampini   hypre_ParCSRMatrix     *hA;
2432df22349SStefano Zampini   hypre_CSRMatrix        *hdiag,*hoffd;
2442df22349SStefano Zampini   MPI_Comm               comm;
2452df22349SStefano Zampini   PetscScalar            *hdd,*hod,*aa,*data;
2467d968826Sstefano_zampini   HYPRE_Int              *col_map_offd,*hdi,*hdj,*hoi,*hoj;
2472df22349SStefano Zampini   PetscInt               *ii,*jj,*iptr,*jptr;
2482df22349SStefano Zampini   PetscInt               cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N;
24958968eb6SStefano Zampini   HYPRE_Int              type;
2502df22349SStefano Zampini   PetscErrorCode         ierr;
2512df22349SStefano Zampini 
2522df22349SStefano Zampini   PetscFunctionBegin;
253a1787963SStefano Zampini   comm = PetscObjectComm((PetscObject)A);
2542df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type));
2552df22349SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
2562df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA));
2572df22349SStefano Zampini   M     = hypre_ParCSRMatrixGlobalNumRows(hA);
2582df22349SStefano Zampini   N     = hypre_ParCSRMatrixGlobalNumCols(hA);
2592df22349SStefano Zampini   str   = hypre_ParCSRMatrixFirstRowIndex(hA);
2602df22349SStefano Zampini   stc   = hypre_ParCSRMatrixFirstColDiag(hA);
2612df22349SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(hA);
2622df22349SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(hA);
2632df22349SStefano Zampini   dr    = hypre_CSRMatrixNumRows(hdiag);
2642df22349SStefano Zampini   dc    = hypre_CSRMatrixNumCols(hdiag);
2652df22349SStefano Zampini   nnz   = hypre_CSRMatrixNumNonzeros(hdiag);
2662df22349SStefano Zampini   hdi   = hypre_CSRMatrixI(hdiag);
2672df22349SStefano Zampini   hdj   = hypre_CSRMatrixJ(hdiag);
2682df22349SStefano Zampini   hdd   = hypre_CSRMatrixData(hdiag);
2692df22349SStefano Zampini   oc    = hypre_CSRMatrixNumCols(hoffd);
2702df22349SStefano Zampini   nnz  += hypre_CSRMatrixNumNonzeros(hoffd);
2712df22349SStefano Zampini   hoi   = hypre_CSRMatrixI(hoffd);
2722df22349SStefano Zampini   hoj   = hypre_CSRMatrixJ(hoffd);
2732df22349SStefano Zampini   hod   = hypre_CSRMatrixData(hoffd);
2742df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
2752df22349SStefano Zampini     PetscInt *aux;
2762df22349SStefano Zampini 
2772df22349SStefano Zampini     /* generate l2g maps for rows and cols */
2782df22349SStefano Zampini     ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr);
2792df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
2802df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2812df22349SStefano Zampini     col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
2822df22349SStefano Zampini     ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr);
2832df22349SStefano Zampini     for (i=0; i<dc; i++) aux[i] = i+stc;
2842df22349SStefano Zampini     for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i];
2852df22349SStefano Zampini     ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
2862df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
2872df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2882df22349SStefano Zampini     /* create MATIS object */
2892df22349SStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
2902df22349SStefano Zampini     ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr);
2912df22349SStefano Zampini     ierr = MatSetType(*B,MATIS);CHKERRQ(ierr);
2922df22349SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr);
2932df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
2942df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
2952df22349SStefano Zampini 
2962df22349SStefano Zampini     /* allocate CSR for local matrix */
2972df22349SStefano Zampini     ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr);
2982df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr);
2992df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr);
3002df22349SStefano Zampini   } else {
3012df22349SStefano Zampini     PetscInt  nr;
3022df22349SStefano Zampini     PetscBool done;
3032df22349SStefano Zampini     ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr);
3042df22349SStefano Zampini     ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr);
3052df22349SStefano 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);
3062df22349SStefano 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);
3072df22349SStefano Zampini     ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr);
3082df22349SStefano Zampini   }
3092df22349SStefano Zampini   /* merge local matrices */
3102df22349SStefano Zampini   ii   = iptr;
3112df22349SStefano Zampini   jj   = jptr;
3122df22349SStefano Zampini   aa   = data;
3132df22349SStefano Zampini   *ii  = *(hdi++) + *(hoi++);
3142df22349SStefano Zampini   for (jd=0,jo=0,cum=0; *ii<nnz; cum++) {
3152df22349SStefano Zampini     PetscScalar *aold = aa;
3162df22349SStefano Zampini     PetscInt    *jold = jj,nc = jd+jo;
3172df22349SStefano Zampini     for (; jd<*hdi; jd++) { *jj++ = *hdj++;      *aa++ = *hdd++; }
3182df22349SStefano Zampini     for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
3192df22349SStefano Zampini     *(++ii) = *(hdi++) + *(hoi++);
3202df22349SStefano Zampini     ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr);
3212df22349SStefano Zampini   }
3222df22349SStefano Zampini   for (; cum<dr; cum++) *(++ii) = nnz;
3232df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
324a033916dSStefano Zampini     Mat_SeqAIJ* a;
325a033916dSStefano Zampini 
3262df22349SStefano Zampini     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr);
3272df22349SStefano Zampini     ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr);
328a033916dSStefano Zampini     /* hack SeqAIJ */
329a033916dSStefano Zampini     a          = (Mat_SeqAIJ*)(lA->data);
330a033916dSStefano Zampini     a->free_a  = PETSC_TRUE;
331a033916dSStefano Zampini     a->free_ij = PETSC_TRUE;
3322df22349SStefano Zampini     ierr = MatDestroy(&lA);CHKERRQ(ierr);
3332df22349SStefano Zampini   }
3342df22349SStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3352df22349SStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3362df22349SStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
3372df22349SStefano Zampini     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
3382df22349SStefano Zampini   }
3392df22349SStefano Zampini   PetscFunctionReturn(0);
3402df22349SStefano Zampini }
3412df22349SStefano Zampini 
3422df22349SStefano Zampini #undef __FUNCT__
34363c07aadSStefano Zampini #define __FUNCT__ "MatConvert_AIJ_HYPRE"
34463c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
34563c07aadSStefano Zampini {
34663c07aadSStefano Zampini   Mat_HYPRE      *hB;
34763c07aadSStefano Zampini   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
34863c07aadSStefano Zampini   PetscErrorCode ierr;
34963c07aadSStefano Zampini 
35063c07aadSStefano Zampini   PetscFunctionBegin;
35163c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
35263c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
35363c07aadSStefano Zampini     /* always destroy the old matrix and create a new memory;
35463c07aadSStefano Zampini        hope this does not churn the memory too much. The problem
35563c07aadSStefano Zampini        is I do not know if it is possible to put the matrix back to
35663c07aadSStefano Zampini        its initial state so that we can directly copy the values
35763c07aadSStefano Zampini        the second time through. */
35863c07aadSStefano Zampini     hB = (Mat_HYPRE*)((*B)->data);
35963c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
36063c07aadSStefano Zampini   } else {
36163c07aadSStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
36263c07aadSStefano Zampini     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
36363c07aadSStefano Zampini     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
36463c07aadSStefano Zampini     hB   = (Mat_HYPRE*)((*B)->data);
36563c07aadSStefano Zampini   }
36663c07aadSStefano Zampini   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
36763c07aadSStefano Zampini   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
3684ec6421dSstefano_zampini   (*B)->preallocated = PETSC_TRUE;
36963c07aadSStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
37063c07aadSStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
37163c07aadSStefano Zampini   PetscFunctionReturn(0);
37263c07aadSStefano Zampini }
37363c07aadSStefano Zampini 
37463c07aadSStefano Zampini #undef __FUNCT__
37563c07aadSStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_AIJ"
376ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
37763c07aadSStefano Zampini {
37863c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
37963c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
38063c07aadSStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
38163c07aadSStefano Zampini   MPI_Comm           comm;
38263c07aadSStefano Zampini   PetscScalar        *da,*oa,*aptr;
38363c07aadSStefano Zampini   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
38463c07aadSStefano Zampini   PetscInt           i,dnnz,onnz,m,n;
38558968eb6SStefano Zampini   HYPRE_Int          type;
38663c07aadSStefano Zampini   PetscMPIInt        size;
38763c07aadSStefano Zampini   PetscErrorCode     ierr;
38863c07aadSStefano Zampini 
38963c07aadSStefano Zampini   PetscFunctionBegin;
39063c07aadSStefano Zampini   comm = PetscObjectComm((PetscObject)A);
39163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
39263c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
39363c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
39463c07aadSStefano Zampini     PetscBool ismpiaij,isseqaij;
39563c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
39663c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
39763c07aadSStefano Zampini     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
39863c07aadSStefano Zampini   }
39963c07aadSStefano Zampini   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
40063c07aadSStefano Zampini 
40163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
40263c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(parcsr);
40363c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(parcsr);
40463c07aadSStefano Zampini   m     = hypre_CSRMatrixNumRows(hdiag);
40563c07aadSStefano Zampini   n     = hypre_CSRMatrixNumCols(hdiag);
40663c07aadSStefano Zampini   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
40763c07aadSStefano Zampini   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
408225daaf8SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
40963c07aadSStefano Zampini     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
41063c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
41163c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
412225daaf8SStefano Zampini   } else if (reuse == MAT_REUSE_MATRIX) {
41363c07aadSStefano Zampini     PetscInt  nr;
41463c07aadSStefano Zampini     PetscBool done;
41563c07aadSStefano Zampini     if (size > 1) {
41663c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
41763c07aadSStefano Zampini 
41863c07aadSStefano Zampini       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
41963c07aadSStefano 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);
42063c07aadSStefano 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);
42163c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
42263c07aadSStefano Zampini     } else {
42363c07aadSStefano Zampini       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
42463c07aadSStefano Zampini       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
42563c07aadSStefano 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);
42663c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
42763c07aadSStefano Zampini     }
428225daaf8SStefano Zampini   } else { /* MAT_INPLACE_MATRIX */
4297d968826Sstefano_zampini     dii = (PetscInt*)hypre_CSRMatrixI(hdiag);
4307d968826Sstefano_zampini     djj = (PetscInt*)hypre_CSRMatrixJ(hdiag);
431225daaf8SStefano Zampini     da  = hypre_CSRMatrixData(hdiag);
43263c07aadSStefano Zampini   }
43363c07aadSStefano Zampini   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
43463c07aadSStefano Zampini   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
43563c07aadSStefano Zampini   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
43663c07aadSStefano Zampini   iptr = djj;
43763c07aadSStefano Zampini   aptr = da;
43863c07aadSStefano Zampini   for (i=0; i<m; i++) {
43963c07aadSStefano Zampini     PetscInt nc = dii[i+1]-dii[i];
44063c07aadSStefano Zampini     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
44163c07aadSStefano Zampini     iptr += nc;
44263c07aadSStefano Zampini     aptr += nc;
44363c07aadSStefano Zampini   }
44463c07aadSStefano Zampini   if (size > 1) {
4457d968826Sstefano_zampini     HYPRE_Int *offdj,*coffd;
44663c07aadSStefano Zampini 
447225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
44863c07aadSStefano Zampini       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
44963c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
45063c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
451225daaf8SStefano Zampini     } else if (reuse == MAT_REUSE_MATRIX) {
45263c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
45363c07aadSStefano Zampini       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
45463c07aadSStefano Zampini       PetscBool  done;
45563c07aadSStefano Zampini 
45663c07aadSStefano Zampini       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
45763c07aadSStefano 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);
45863c07aadSStefano 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);
45963c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
460225daaf8SStefano Zampini     } else { /* MAT_INPLACE_MATRIX */
4617d968826Sstefano_zampini       oii = (PetscInt*)hypre_CSRMatrixI(hoffd);
4627d968826Sstefano_zampini       ojj = (PetscInt*)hypre_CSRMatrixJ(hoffd);
463225daaf8SStefano Zampini       oa  = hypre_CSRMatrixData(hoffd);
46463c07aadSStefano Zampini     }
46563c07aadSStefano Zampini     ierr  = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
46663c07aadSStefano Zampini     offdj = hypre_CSRMatrixJ(hoffd);
46763c07aadSStefano Zampini     coffd = hypre_ParCSRMatrixColMapOffd(parcsr);
46863c07aadSStefano Zampini     for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]];
46963c07aadSStefano Zampini     ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr);
47063c07aadSStefano Zampini     iptr = ojj;
47163c07aadSStefano Zampini     aptr = oa;
47263c07aadSStefano Zampini     for (i=0; i<m; i++) {
47363c07aadSStefano Zampini        PetscInt nc = oii[i+1]-oii[i];
47463c07aadSStefano Zampini        ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
47563c07aadSStefano Zampini        iptr += nc;
47663c07aadSStefano Zampini        aptr += nc;
47763c07aadSStefano Zampini     }
478225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
47963c07aadSStefano Zampini       Mat_MPIAIJ *b;
48063c07aadSStefano Zampini       Mat_SeqAIJ *d,*o;
481225daaf8SStefano Zampini 
48241073d71Sstefano_zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,B);CHKERRQ(ierr);
48363c07aadSStefano Zampini       /* hack MPIAIJ */
48463c07aadSStefano Zampini       b          = (Mat_MPIAIJ*)((*B)->data);
48563c07aadSStefano Zampini       d          = (Mat_SeqAIJ*)b->A->data;
48663c07aadSStefano Zampini       o          = (Mat_SeqAIJ*)b->B->data;
48763c07aadSStefano Zampini       d->free_a  = PETSC_TRUE;
48863c07aadSStefano Zampini       d->free_ij = PETSC_TRUE;
48963c07aadSStefano Zampini       o->free_a  = PETSC_TRUE;
49063c07aadSStefano Zampini       o->free_ij = PETSC_TRUE;
491225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
492225daaf8SStefano Zampini       Mat T;
49341073d71Sstefano_zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,&T);CHKERRQ(ierr);
494225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
495225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
496225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
497225daaf8SStefano Zampini       hypre_CSRMatrixI(hoffd)    = NULL;
498225daaf8SStefano Zampini       hypre_CSRMatrixJ(hoffd)    = NULL;
499225daaf8SStefano Zampini       hypre_CSRMatrixData(hoffd) = NULL;
500225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
50163c07aadSStefano Zampini     }
502225daaf8SStefano Zampini   } else {
503225daaf8SStefano Zampini     oii  = NULL;
504225daaf8SStefano Zampini     ojj  = NULL;
505225daaf8SStefano Zampini     oa   = NULL;
506225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
50763c07aadSStefano Zampini       Mat_SeqAIJ* b;
50863c07aadSStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
50963c07aadSStefano Zampini       /* hack SeqAIJ */
51063c07aadSStefano Zampini       b          = (Mat_SeqAIJ*)((*B)->data);
51163c07aadSStefano Zampini       b->free_a  = PETSC_TRUE;
51263c07aadSStefano Zampini       b->free_ij = PETSC_TRUE;
513225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
514225daaf8SStefano Zampini       Mat T;
515225daaf8SStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,&T);CHKERRQ(ierr);
516225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
517225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
518225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
519225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
52063c07aadSStefano Zampini     }
521225daaf8SStefano Zampini   }
522225daaf8SStefano Zampini 
523225daaf8SStefano Zampini   /* we have to use hypre_Tfree to free the arrays */
52463c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
525225daaf8SStefano Zampini     void *ptrs[6] = {dii,djj,da,oii,ojj,oa};
526225daaf8SStefano Zampini     const char *names[6] = {"_hypre_csr_dii",
527225daaf8SStefano Zampini                             "_hypre_csr_djj",
528225daaf8SStefano Zampini                             "_hypre_csr_da",
529225daaf8SStefano Zampini                             "_hypre_csr_oii",
530225daaf8SStefano Zampini                             "_hypre_csr_ojj",
531225daaf8SStefano Zampini                             "_hypre_csr_oa"};
532225daaf8SStefano Zampini     for (i=0; i<6; i++) {
533225daaf8SStefano Zampini       PetscContainer c;
534225daaf8SStefano Zampini 
535225daaf8SStefano Zampini       ierr = PetscContainerCreate(comm,&c);CHKERRQ(ierr);
536225daaf8SStefano Zampini       ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr);
537225daaf8SStefano Zampini       ierr = PetscContainerSetUserDestroy(c,hypre_array_destroy);CHKERRQ(ierr);
538225daaf8SStefano Zampini       ierr = PetscObjectCompose((PetscObject)(*B),names[i],(PetscObject)c);CHKERRQ(ierr);
539225daaf8SStefano Zampini       ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
540225daaf8SStefano Zampini     }
54163c07aadSStefano Zampini   }
54263c07aadSStefano Zampini   PetscFunctionReturn(0);
54363c07aadSStefano Zampini }
54463c07aadSStefano Zampini 
54563c07aadSStefano Zampini #undef __FUNCT__
546613e5ff0Sstefano_zampini #define __FUNCT__ "MatAIJGetParCSR_Private"
547613e5ff0Sstefano_zampini static PetscErrorCode MatAIJGetParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
548c1a070e6SStefano Zampini {
549613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *tA;
550c1a070e6SStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
551c1a070e6SStefano Zampini   Mat_SeqAIJ         *diag,*offd;
552c1a070e6SStefano Zampini   PetscInt           *garray,noffd,dnnz,onnz,*row_starts,*col_starts;
553c1a070e6SStefano Zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
554613e5ff0Sstefano_zampini   PetscBool          ismpiaij,isseqaij;
555c1a070e6SStefano Zampini   PetscErrorCode     ierr;
556c1a070e6SStefano Zampini 
557c1a070e6SStefano Zampini   PetscFunctionBegin;
558c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
559c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
560c1a070e6SStefano Zampini   if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type);
561c1a070e6SStefano Zampini   if (ismpiaij) {
562c1a070e6SStefano Zampini     Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data);
563c1a070e6SStefano Zampini 
564c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)a->A->data;
565c1a070e6SStefano Zampini     offd   = (Mat_SeqAIJ*)a->B->data;
566c1a070e6SStefano Zampini     garray = a->garray;
567c1a070e6SStefano Zampini     noffd  = a->B->cmap->N;
568c1a070e6SStefano Zampini     dnnz   = diag->nz;
569c1a070e6SStefano Zampini     onnz   = offd->nz;
570c1a070e6SStefano Zampini   } else {
571c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)A->data;
572c1a070e6SStefano Zampini     offd   = NULL;
573c1a070e6SStefano Zampini     garray = NULL;
574c1a070e6SStefano Zampini     noffd  = 0;
575c1a070e6SStefano Zampini     dnnz   = diag->nz;
576c1a070e6SStefano Zampini     onnz   = 0;
577c1a070e6SStefano Zampini   }
578225daaf8SStefano Zampini 
579c1a070e6SStefano Zampini   /* create a temporary ParCSR */
580c1a070e6SStefano Zampini   if (HYPRE_AssumedPartitionCheck()) {
581c1a070e6SStefano Zampini     PetscMPIInt myid;
582c1a070e6SStefano Zampini 
583c1a070e6SStefano Zampini     ierr       = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
584c1a070e6SStefano Zampini     row_starts = A->rmap->range + myid;
585c1a070e6SStefano Zampini     col_starts = A->cmap->range + myid;
586c1a070e6SStefano Zampini   } else {
587c1a070e6SStefano Zampini     row_starts = A->rmap->range;
588c1a070e6SStefano Zampini     col_starts = A->cmap->range;
589c1a070e6SStefano Zampini   }
5907d968826Sstefano_zampini   tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,(HYPRE_Int*)row_starts,(HYPRE_Int*)col_starts,noffd,dnnz,onnz);
591c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(tA,0);
592c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(tA,0);
593c1a070e6SStefano Zampini 
594225daaf8SStefano Zampini   /* set diagonal part */
595c1a070e6SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(tA);
5967d968826Sstefano_zampini   hypre_CSRMatrixI(hdiag)           = (HYPRE_Int*)diag->i;
5977d968826Sstefano_zampini   hypre_CSRMatrixJ(hdiag)           = (HYPRE_Int*)diag->j;
598c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = diag->a;
599c1a070e6SStefano Zampini   hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz;
600c1a070e6SStefano Zampini   hypre_CSRMatrixSetRownnz(hdiag);
601c1a070e6SStefano Zampini   hypre_CSRMatrixSetDataOwner(hdiag,0);
602c1a070e6SStefano Zampini 
603225daaf8SStefano Zampini   /* set offdiagonal part */
604c1a070e6SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(tA);
605c1a070e6SStefano Zampini   if (offd) {
6067d968826Sstefano_zampini     hypre_CSRMatrixI(hoffd)           = (HYPRE_Int*)offd->i;
6077d968826Sstefano_zampini     hypre_CSRMatrixJ(hoffd)           = (HYPRE_Int*)offd->j;
608c1a070e6SStefano Zampini     hypre_CSRMatrixData(hoffd)        = offd->a;
609c1a070e6SStefano Zampini     hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz;
610c1a070e6SStefano Zampini     hypre_CSRMatrixSetRownnz(hoffd);
611c1a070e6SStefano Zampini     hypre_CSRMatrixSetDataOwner(hoffd,0);
612c1a070e6SStefano Zampini     hypre_ParCSRMatrixSetNumNonzeros(tA);
6137d968826Sstefano_zampini     hypre_ParCSRMatrixColMapOffd(tA) = (HYPRE_Int*)garray;
614c1a070e6SStefano Zampini   }
615613e5ff0Sstefano_zampini   *hA = tA;
616613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
617613e5ff0Sstefano_zampini }
618c1a070e6SStefano Zampini 
619613e5ff0Sstefano_zampini #undef __FUNCT__
620613e5ff0Sstefano_zampini #define __FUNCT__ "MatAIJRestoreParCSR_Private"
621613e5ff0Sstefano_zampini static PetscErrorCode MatAIJRestoreParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
622613e5ff0Sstefano_zampini {
623613e5ff0Sstefano_zampini   hypre_CSRMatrix    *hdiag,*hoffd;
624c1a070e6SStefano Zampini 
625613e5ff0Sstefano_zampini   PetscFunctionBegin;
626613e5ff0Sstefano_zampini   hdiag = hypre_ParCSRMatrixDiag(*hA);
627613e5ff0Sstefano_zampini   hoffd = hypre_ParCSRMatrixOffd(*hA);
628c1a070e6SStefano Zampini   /* set pointers to NULL before destroying tA */
629c1a070e6SStefano Zampini   hypre_CSRMatrixI(hdiag)           = NULL;
630c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hdiag)           = NULL;
631c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = NULL;
632c1a070e6SStefano Zampini   hypre_CSRMatrixI(hoffd)           = NULL;
633c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hoffd)           = NULL;
634c1a070e6SStefano Zampini   hypre_CSRMatrixData(hoffd)        = NULL;
635613e5ff0Sstefano_zampini   hypre_ParCSRMatrixColMapOffd(*hA) = NULL;
636613e5ff0Sstefano_zampini   hypre_ParCSRMatrixDestroy(*hA);
637613e5ff0Sstefano_zampini   *hA = NULL;
638613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
639613e5ff0Sstefano_zampini }
640613e5ff0Sstefano_zampini 
641613e5ff0Sstefano_zampini /* calls RAP from BoomerAMG:
642*3dad0653Sstefano_zampini    the resulting ParCSR will not own the column and row starts
643*3dad0653Sstefano_zampini    It looks like we don't need to have the diagonal entries
644*3dad0653Sstefano_zampini    ordered first in the rows of the diagonal part
645*3dad0653Sstefano_zampini    for boomerAMGBuildCoarseOperator to work */
646613e5ff0Sstefano_zampini #undef __FUNCT__
647613e5ff0Sstefano_zampini #define __FUNCT__ "MatHYPRE_ParCSR_RAP"
648613e5ff0Sstefano_zampini static PetscErrorCode MatHYPRE_ParCSR_RAP(hypre_ParCSRMatrix *hR, hypre_ParCSRMatrix *hA,
649613e5ff0Sstefano_zampini                                           hypre_ParCSRMatrix *hP, hypre_ParCSRMatrix **hRAP)
650613e5ff0Sstefano_zampini {
651613e5ff0Sstefano_zampini   PetscErrorCode ierr;
652613e5ff0Sstefano_zampini   HYPRE_Int      P_owns_col_starts,R_owns_row_starts;
653613e5ff0Sstefano_zampini 
654613e5ff0Sstefano_zampini   PetscFunctionBegin;
655613e5ff0Sstefano_zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(hP);
656613e5ff0Sstefano_zampini   R_owns_row_starts = hypre_ParCSRMatrixOwnsRowStarts(hR);
657613e5ff0Sstefano_zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(hR,hA,hP,hRAP));
658613e5ff0Sstefano_zampini   PetscStackCallStandard(hypre_ParCSRMatrixSetNumNonzeros,(*hRAP));
659613e5ff0Sstefano_zampini   /* hypre_BoomerAMGBuildCoarseOperator steals the col_starts from P and the row_starts from R */
660613e5ff0Sstefano_zampini   hypre_ParCSRMatrixSetRowStartsOwner(*hRAP,0);
661613e5ff0Sstefano_zampini   hypre_ParCSRMatrixSetColStartsOwner(*hRAP,0);
662613e5ff0Sstefano_zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(hP,1);
663613e5ff0Sstefano_zampini   if (R_owns_row_starts) hypre_ParCSRMatrixSetRowStartsOwner(hR,1);
664613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
665613e5ff0Sstefano_zampini }
666613e5ff0Sstefano_zampini 
667613e5ff0Sstefano_zampini #undef __FUNCT__
6686f231fbdSstefano_zampini #define __FUNCT__ "MatPtAPNumeric_AIJ_AIJ_wHYPRE"
6696f231fbdSstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat P,Mat C)
670613e5ff0Sstefano_zampini {
6716f231fbdSstefano_zampini   Mat                B;
672613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *hA,*hP,*hPtAP;
673613e5ff0Sstefano_zampini   PetscErrorCode     ierr;
674613e5ff0Sstefano_zampini 
675613e5ff0Sstefano_zampini   PetscFunctionBegin;
676613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
677613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(P,&hP);CHKERRQ(ierr);
678613e5ff0Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(hP,hA,hP,&hPtAP);CHKERRQ(ierr);
6796f231fbdSstefano_zampini   ierr = MatCreateFromParCSR(hPtAP,MATAIJ,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
6806f231fbdSstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
681613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
682613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(P,&hP);CHKERRQ(ierr);
6836f231fbdSstefano_zampini   PetscFunctionReturn(0);
6846f231fbdSstefano_zampini }
6856f231fbdSstefano_zampini 
6866f231fbdSstefano_zampini #undef __FUNCT__
6876f231fbdSstefano_zampini #define __FUNCT__ "MatPtAPSymbolic_AIJ_AIJ_wHYPRE"
6886f231fbdSstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat P,PetscReal fill,Mat *C)
6896f231fbdSstefano_zampini {
6906f231fbdSstefano_zampini   PetscErrorCode     ierr;
6916f231fbdSstefano_zampini   PetscFunctionBegin;
6926f231fbdSstefano_zampini   ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
6936f231fbdSstefano_zampini   ierr                   = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
6946f231fbdSstefano_zampini   (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_AIJ_wHYPRE;
695613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
696613e5ff0Sstefano_zampini }
697613e5ff0Sstefano_zampini 
698613e5ff0Sstefano_zampini #undef __FUNCT__
6994cc28894Sstefano_zampini #define __FUNCT__ "MatPtAPNumeric_AIJ_HYPRE"
7004cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_HYPRE(Mat A,Mat P,Mat C)
701613e5ff0Sstefano_zampini {
7024cc28894Sstefano_zampini   Mat                B;
7034cc28894Sstefano_zampini   Mat_HYPRE          *hP;
704613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *hA = NULL,*Pparcsr,*ptapparcsr;
705613e5ff0Sstefano_zampini   HYPRE_Int          type;
706613e5ff0Sstefano_zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
7074cc28894Sstefano_zampini   PetscBool          ishypre;
708613e5ff0Sstefano_zampini   PetscErrorCode     ierr;
709613e5ff0Sstefano_zampini 
710613e5ff0Sstefano_zampini   PetscFunctionBegin;
7114cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
7124cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(comm,PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
7134cc28894Sstefano_zampini   hP = (Mat_HYPRE*)P->data;
714613e5ff0Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
715613e5ff0Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
716613e5ff0Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
717613e5ff0Sstefano_zampini 
718613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
719613e5ff0Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,hA,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
720613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
721225daaf8SStefano Zampini 
7224cc28894Sstefano_zampini   /* create temporary matrix and merge to C */
7234cc28894Sstefano_zampini   ierr = MatCreateFromParCSR(ptapparcsr,((PetscObject)C)->type_name,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
7244cc28894Sstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
7254cc28894Sstefano_zampini   PetscFunctionReturn(0);
7264cc28894Sstefano_zampini }
7274cc28894Sstefano_zampini 
7284cc28894Sstefano_zampini #undef __FUNCT__
7294cc28894Sstefano_zampini #define __FUNCT__ "MatPtAP_AIJ_HYPRE"
7304cc28894Sstefano_zampini static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
7314cc28894Sstefano_zampini {
7324cc28894Sstefano_zampini   PetscErrorCode ierr;
7334cc28894Sstefano_zampini 
7344cc28894Sstefano_zampini   PetscFunctionBegin;
7354cc28894Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
7364cc28894Sstefano_zampini     const char *deft = MATAIJ;
7374cc28894Sstefano_zampini     char       type[256];
7384cc28894Sstefano_zampini     PetscBool  flg;
7394cc28894Sstefano_zampini 
7404cc28894Sstefano_zampini     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
7414cc28894Sstefano_zampini     ierr = PetscOptionsFList("-matptap_hypre_outtype","Matrix type",NULL,MatList,deft,type,256,&flg);CHKERRQ(ierr);
7424cc28894Sstefano_zampini     ierr = PetscOptionsEnd();CHKERRQ(ierr);
7434cc28894Sstefano_zampini     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7444cc28894Sstefano_zampini     ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
7454cc28894Sstefano_zampini     if (flg) {
7464cc28894Sstefano_zampini       ierr = MatSetType(*C,type);CHKERRQ(ierr);
7474cc28894Sstefano_zampini     } else {
7484cc28894Sstefano_zampini       ierr = MatSetType(*C,deft);CHKERRQ(ierr);
7494cc28894Sstefano_zampini     }
7504cc28894Sstefano_zampini     (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_HYPRE;
7514cc28894Sstefano_zampini     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7524cc28894Sstefano_zampini   }
7534cc28894Sstefano_zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
7544cc28894Sstefano_zampini   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
755613e5ff0Sstefano_zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
756c1a070e6SStefano Zampini   PetscFunctionReturn(0);
757c1a070e6SStefano Zampini }
758c1a070e6SStefano Zampini 
759c1a070e6SStefano Zampini #undef __FUNCT__
7604cc28894Sstefano_zampini #define __FUNCT__ "MatPtAPNumeric_HYPRE_HYPRE"
7614cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_HYPRE_HYPRE(Mat A,Mat P,Mat C)
7624cc28894Sstefano_zampini {
7634cc28894Sstefano_zampini   Mat                B;
7644cc28894Sstefano_zampini   hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr;
7654cc28894Sstefano_zampini   Mat_HYPRE          *hA,*hP;
7664cc28894Sstefano_zampini   PetscBool          ishypre;
7674cc28894Sstefano_zampini   HYPRE_Int          type;
7684cc28894Sstefano_zampini   PetscErrorCode     ierr;
7694cc28894Sstefano_zampini 
7704cc28894Sstefano_zampini   PetscFunctionBegin;
7714cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
7724cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)P),PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
7734cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr);
7744cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE);
7754cc28894Sstefano_zampini   hA = (Mat_HYPRE*)A->data;
7764cc28894Sstefano_zampini   hP = (Mat_HYPRE*)P->data;
7774cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
7784cc28894Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
7794cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
7804cc28894Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
7814cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
7824cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
7834cc28894Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
7844cc28894Sstefano_zampini   ierr = MatCreateFromParCSR(ptapparcsr,MATHYPRE,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
7854cc28894Sstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
7864cc28894Sstefano_zampini   PetscFunctionReturn(0);
7874cc28894Sstefano_zampini }
7884cc28894Sstefano_zampini 
7894cc28894Sstefano_zampini #undef __FUNCT__
790cd8bc7baSStefano Zampini #define __FUNCT__ "MatPtAP_HYPRE_HYPRE"
791cd8bc7baSStefano Zampini static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
792cd8bc7baSStefano Zampini {
793cd8bc7baSStefano Zampini   PetscErrorCode ierr;
794cd8bc7baSStefano Zampini 
795cd8bc7baSStefano Zampini   PetscFunctionBegin;
7964cc28894Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
7974cc28894Sstefano_zampini     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7984cc28894Sstefano_zampini     ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
7994cc28894Sstefano_zampini     ierr                   = MatSetType(*C,MATHYPRE);CHKERRQ(ierr);
8004cc28894Sstefano_zampini     (*C)->ops->ptapnumeric = MatPtAPNumeric_HYPRE_HYPRE;
8014cc28894Sstefano_zampini     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
8024cc28894Sstefano_zampini   }
803613e5ff0Sstefano_zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
8044cc28894Sstefano_zampini   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
805613e5ff0Sstefano_zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
806cd8bc7baSStefano Zampini   PetscFunctionReturn(0);
807cd8bc7baSStefano Zampini }
808cd8bc7baSStefano Zampini 
809d501dc42Sstefano_zampini /* calls hypre_ParMatmul
810d501dc42Sstefano_zampini    hypre_ParMatMul uses hypre_ParMatrixCreate with the communicator of hA
811*3dad0653Sstefano_zampini    hypre_ParMatrixCreate does not duplicate the communicator
812*3dad0653Sstefano_zampini    It looks like we don't need to have the diagonal entries
813*3dad0653Sstefano_zampini    ordered first in the rows of the diagonal part
814*3dad0653Sstefano_zampini    for boomerAMGBuildCoarseOperator to work */
815d501dc42Sstefano_zampini #undef __FUNCT__
816d501dc42Sstefano_zampini #define __FUNCT__ "MatHYPRE_ParCSR_MatMatMult"
817d501dc42Sstefano_zampini static PetscErrorCode MatHYPRE_ParCSR_MatMatMult(hypre_ParCSRMatrix *hA, hypre_ParCSRMatrix *hB, hypre_ParCSRMatrix **hAB)
818d501dc42Sstefano_zampini {
819d501dc42Sstefano_zampini   PetscFunctionBegin;
820d501dc42Sstefano_zampini   PetscStackPush("hypre_ParMatmul");
821d501dc42Sstefano_zampini   *hAB = hypre_ParMatmul(hA,hB);
822d501dc42Sstefano_zampini   PetscStackPop;
823d501dc42Sstefano_zampini   PetscFunctionReturn(0);
824d501dc42Sstefano_zampini }
825d501dc42Sstefano_zampini 
826cd8bc7baSStefano Zampini #undef __FUNCT__
8275e5acdf2Sstefano_zampini #define __FUNCT__ "MatMatMultNumeric_AIJ_AIJ_wHYPRE"
8285e5acdf2Sstefano_zampini static PetscErrorCode MatMatMultNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C)
8295e5acdf2Sstefano_zampini {
8305e5acdf2Sstefano_zampini   Mat                D;
831d501dc42Sstefano_zampini   hypre_ParCSRMatrix *hA,*hB,*hAB = NULL;
8325e5acdf2Sstefano_zampini   PetscErrorCode     ierr;
8335e5acdf2Sstefano_zampini 
8345e5acdf2Sstefano_zampini   PetscFunctionBegin;
8355e5acdf2Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
8365e5acdf2Sstefano_zampini   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
837d501dc42Sstefano_zampini   ierr = MatHYPRE_ParCSR_MatMatMult(hA,hB,&hAB);CHKERRQ(ierr);
8385e5acdf2Sstefano_zampini   ierr = MatCreateFromParCSR(hAB,MATAIJ,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
8395e5acdf2Sstefano_zampini   ierr = MatHeaderMerge(C,&D);CHKERRQ(ierr);
8405e5acdf2Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
8415e5acdf2Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
8425e5acdf2Sstefano_zampini   PetscFunctionReturn(0);
8435e5acdf2Sstefano_zampini }
8445e5acdf2Sstefano_zampini 
8455e5acdf2Sstefano_zampini #undef __FUNCT__
8465e5acdf2Sstefano_zampini #define __FUNCT__ "MatMatMultSymbolic_AIJ_AIJ_wHYPRE"
8475e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat B,PetscReal fill,Mat *C)
8485e5acdf2Sstefano_zampini {
8495e5acdf2Sstefano_zampini   PetscErrorCode ierr;
85020e1dc0dSstefano_zampini 
8515e5acdf2Sstefano_zampini   PetscFunctionBegin;
8525e5acdf2Sstefano_zampini   ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
8535e5acdf2Sstefano_zampini   ierr                      = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
8545e5acdf2Sstefano_zampini   (*C)->ops->matmultnumeric = MatMatMultNumeric_AIJ_AIJ_wHYPRE;
8555e5acdf2Sstefano_zampini   PetscFunctionReturn(0);
8565e5acdf2Sstefano_zampini }
8575e5acdf2Sstefano_zampini 
8585e5acdf2Sstefano_zampini #undef __FUNCT__
859d501dc42Sstefano_zampini #define __FUNCT__ "MatMatMultNumeric_HYPRE_HYPRE"
860d501dc42Sstefano_zampini static PetscErrorCode MatMatMultNumeric_HYPRE_HYPRE(Mat A,Mat B,Mat C)
861d501dc42Sstefano_zampini {
862d501dc42Sstefano_zampini   Mat                D;
863d501dc42Sstefano_zampini   hypre_ParCSRMatrix *Aparcsr,*Bparcsr,*ABparcsr = NULL;
864d501dc42Sstefano_zampini   Mat_HYPRE          *hA,*hB;
865d501dc42Sstefano_zampini   PetscBool          ishypre;
866d501dc42Sstefano_zampini   HYPRE_Int          type;
867d501dc42Sstefano_zampini   PetscErrorCode     ierr;
868d501dc42Sstefano_zampini 
869d501dc42Sstefano_zampini   PetscFunctionBegin;
870d501dc42Sstefano_zampini   ierr = PetscPrintf(PetscObjectComm((PetscObject)A),"USING %s\n",__FUNCT__);CHKERRQ(ierr);
871d501dc42Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)B,MATHYPRE,&ishypre);CHKERRQ(ierr);
872d501dc42Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_USER,"B should be of type %s",MATHYPRE);
873d501dc42Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr);
874d501dc42Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE);
875d501dc42Sstefano_zampini   hA = (Mat_HYPRE*)A->data;
876d501dc42Sstefano_zampini   hB = (Mat_HYPRE*)B->data;
877d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
878d501dc42Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
879d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hB->ij,&type));
880d501dc42Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
881d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
882d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hB->ij,(void**)&Bparcsr));
883d501dc42Sstefano_zampini   ierr = MatHYPRE_ParCSR_MatMatMult(Aparcsr,Bparcsr,&ABparcsr);CHKERRQ(ierr);
884d501dc42Sstefano_zampini   ierr = MatCreateFromParCSR(ABparcsr,MATHYPRE,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
885d501dc42Sstefano_zampini   /* need to use HeaderReplace because HeaderMerge messes up with the communicator */
886d501dc42Sstefano_zampini   ierr = MatHeaderReplace(C,&D);CHKERRQ(ierr);
887d501dc42Sstefano_zampini   C->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE;
888d501dc42Sstefano_zampini   PetscFunctionReturn(0);
889d501dc42Sstefano_zampini }
890d501dc42Sstefano_zampini 
891d501dc42Sstefano_zampini #undef __FUNCT__
892d501dc42Sstefano_zampini #define __FUNCT__ "MatMatMult_HYPRE_HYPRE"
893d501dc42Sstefano_zampini static PetscErrorCode MatMatMult_HYPRE_HYPRE(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
894d501dc42Sstefano_zampini {
895d501dc42Sstefano_zampini   PetscErrorCode ierr;
896d501dc42Sstefano_zampini 
897d501dc42Sstefano_zampini   PetscFunctionBegin;
898d501dc42Sstefano_zampini   ierr = PetscPrintf(PetscObjectComm((PetscObject)A),"USING %s\n",__FUNCT__);CHKERRQ(ierr);
899d501dc42Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
900d501dc42Sstefano_zampini     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
901d501dc42Sstefano_zampini     ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
902d501dc42Sstefano_zampini     ierr                      = MatSetType(*C,MATHYPRE);CHKERRQ(ierr);
903d501dc42Sstefano_zampini     (*C)->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE;
904d501dc42Sstefano_zampini     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
905d501dc42Sstefano_zampini   }
906d501dc42Sstefano_zampini   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
907d501dc42Sstefano_zampini   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
908d501dc42Sstefano_zampini   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
909d501dc42Sstefano_zampini   PetscFunctionReturn(0);
910d501dc42Sstefano_zampini }
911d501dc42Sstefano_zampini 
912d501dc42Sstefano_zampini 
913d501dc42Sstefano_zampini #undef __FUNCT__
914*3dad0653Sstefano_zampini #define __FUNCT__ "MatTransposeMatMatMultNumeric_AIJ_AIJ_AIJ_wHYPRE"
915*3dad0653Sstefano_zampini PETSC_INTERN PetscErrorCode MatTransposeMatMatMultNumeric_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,Mat D)
91620e1dc0dSstefano_zampini {
91720e1dc0dSstefano_zampini   Mat                E;
91820e1dc0dSstefano_zampini   hypre_ParCSRMatrix *hA,*hB,*hC,*hABC;
91920e1dc0dSstefano_zampini   PetscErrorCode     ierr;
92020e1dc0dSstefano_zampini 
92120e1dc0dSstefano_zampini   PetscFunctionBegin;
92220e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
92320e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
92420e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(C,&hC);CHKERRQ(ierr);
92520e1dc0dSstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(hA,hB,hC,&hABC);CHKERRQ(ierr);
92620e1dc0dSstefano_zampini   ierr = MatCreateFromParCSR(hABC,MATAIJ,PETSC_OWN_POINTER,&E);CHKERRQ(ierr);
92720e1dc0dSstefano_zampini   ierr = MatHeaderMerge(D,&E);CHKERRQ(ierr);
92820e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
92920e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
93020e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(C,&hC);CHKERRQ(ierr);
93120e1dc0dSstefano_zampini   PetscFunctionReturn(0);
93220e1dc0dSstefano_zampini }
93320e1dc0dSstefano_zampini 
93420e1dc0dSstefano_zampini #undef __FUNCT__
935*3dad0653Sstefano_zampini #define __FUNCT__ "MatTransposeMatMatMultSymbolic_AIJ_AIJ_AIJ_wHYPRE"
936*3dad0653Sstefano_zampini PETSC_INTERN PetscErrorCode MatTransposeMatMatMultSymbolic_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,PetscReal fill,Mat *D)
93720e1dc0dSstefano_zampini {
93820e1dc0dSstefano_zampini   PetscErrorCode ierr;
93920e1dc0dSstefano_zampini 
94020e1dc0dSstefano_zampini   PetscFunctionBegin;
94120e1dc0dSstefano_zampini   ierr = MatCreate(PetscObjectComm((PetscObject)A),D);CHKERRQ(ierr);
94220e1dc0dSstefano_zampini   ierr = MatSetType(*D,MATAIJ);CHKERRQ(ierr);
94320e1dc0dSstefano_zampini   PetscFunctionReturn(0);
94420e1dc0dSstefano_zampini }
94520e1dc0dSstefano_zampini 
94620e1dc0dSstefano_zampini #undef __FUNCT__
94763c07aadSStefano Zampini #define __FUNCT__ "MatMultTranspose_HYPRE"
948ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
94963c07aadSStefano Zampini {
95063c07aadSStefano Zampini   PetscErrorCode ierr;
95163c07aadSStefano Zampini 
95263c07aadSStefano Zampini   PetscFunctionBegin;
95363c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
95463c07aadSStefano Zampini   PetscFunctionReturn(0);
95563c07aadSStefano Zampini }
95663c07aadSStefano Zampini 
95763c07aadSStefano Zampini #undef __FUNCT__
95863c07aadSStefano Zampini #define __FUNCT__ "MatMult_HYPRE"
959ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
96063c07aadSStefano Zampini {
96163c07aadSStefano Zampini   PetscErrorCode ierr;
96263c07aadSStefano Zampini 
96363c07aadSStefano Zampini   PetscFunctionBegin;
96463c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
96563c07aadSStefano Zampini   PetscFunctionReturn(0);
96663c07aadSStefano Zampini }
96763c07aadSStefano Zampini 
96863c07aadSStefano Zampini #undef __FUNCT__
96963c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_MultKernel_Private"
970ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
97163c07aadSStefano Zampini {
97263c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
97363c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
97463c07aadSStefano Zampini   hypre_ParVector    *hx,*hy;
97563c07aadSStefano Zampini   PetscScalar        *ax,*ay,*sax,*say;
97663c07aadSStefano Zampini   PetscErrorCode     ierr;
97763c07aadSStefano Zampini 
97863c07aadSStefano Zampini   PetscFunctionBegin;
97963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
98063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
98163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
98263c07aadSStefano Zampini   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
98363c07aadSStefano Zampini   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
98463c07aadSStefano Zampini   if (trans) {
98558968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ay,say);
98658968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax);
98763c07aadSStefano Zampini     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
98858968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,say,ay);
98958968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax);
99063c07aadSStefano Zampini   } else {
99158968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax);
99258968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ay,say);
99363c07aadSStefano Zampini     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
99458968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax);
99558968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,say,ay);
99663c07aadSStefano Zampini   }
99763c07aadSStefano Zampini   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
99863c07aadSStefano Zampini   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
99963c07aadSStefano Zampini   PetscFunctionReturn(0);
100063c07aadSStefano Zampini }
100163c07aadSStefano Zampini 
100263c07aadSStefano Zampini #undef __FUNCT__
100363c07aadSStefano Zampini #define __FUNCT__ "MatDestroy_HYPRE"
1004ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A)
100563c07aadSStefano Zampini {
100663c07aadSStefano Zampini   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
100763c07aadSStefano Zampini   PetscErrorCode ierr;
100863c07aadSStefano Zampini 
100963c07aadSStefano Zampini   PetscFunctionBegin;
101063c07aadSStefano Zampini   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
101163c07aadSStefano Zampini   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
1012978814f1SStefano Zampini   if (hA->ij) {
1013978814f1SStefano Zampini     if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL;
1014978814f1SStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
1015978814f1SStefano Zampini   }
101663c07aadSStefano Zampini   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); }
101763c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
10182df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr);
1019c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr);
1020c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr);
1021d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPRESetPreallocation_C",NULL);CHKERRQ(ierr);
1022dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPREGetParCSR_C",NULL);CHKERRQ(ierr);
102363c07aadSStefano Zampini   ierr = PetscFree(A->data);CHKERRQ(ierr);
102463c07aadSStefano Zampini   PetscFunctionReturn(0);
102563c07aadSStefano Zampini }
102663c07aadSStefano Zampini 
102763c07aadSStefano Zampini #undef __FUNCT__
102863c07aadSStefano Zampini #define __FUNCT__ "MatSetUp_HYPRE"
1029ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A)
103063c07aadSStefano Zampini {
10314ec6421dSstefano_zampini   PetscErrorCode ierr;
10324ec6421dSstefano_zampini 
10334ec6421dSstefano_zampini   PetscFunctionBegin;
10344ec6421dSstefano_zampini   ierr = MatHYPRESetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr);
10354ec6421dSstefano_zampini   PetscFunctionReturn(0);
10364ec6421dSstefano_zampini }
10374ec6421dSstefano_zampini 
10384ec6421dSstefano_zampini #undef __FUNCT__
10394ec6421dSstefano_zampini #define __FUNCT__ "MatAssemblyEnd_HYPRE"
10404ec6421dSstefano_zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
10414ec6421dSstefano_zampini {
104263c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
104363c07aadSStefano Zampini   Vec                x,b;
104463c07aadSStefano Zampini   PetscErrorCode     ierr;
104563c07aadSStefano Zampini 
104663c07aadSStefano Zampini   PetscFunctionBegin;
10474ec6421dSstefano_zampini   if (mode == MAT_FLUSH_ASSEMBLY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_FLUSH_ASSEMBLY currently not supported with MATHYPRE");
10484ec6421dSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
104963c07aadSStefano Zampini   if (hA->x) PetscFunctionReturn(0);
105063c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
105163c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
105263c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
105363c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
105463c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
105563c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
105663c07aadSStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
105763c07aadSStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
105863c07aadSStefano Zampini   PetscFunctionReturn(0);
105963c07aadSStefano Zampini }
106063c07aadSStefano Zampini 
1061d975228cSstefano_zampini #define MATHYPRE_SCRATCH 2048
1062d975228cSstefano_zampini 
1063d975228cSstefano_zampini #undef __FUNCT__
1064d975228cSstefano_zampini #define __FUNCT__ "MatSetValues_HYPRE"
1065d975228cSstefano_zampini PetscErrorCode MatSetValues_HYPRE(Mat A, PetscInt nr, const PetscInt rows[], PetscInt nc, const PetscInt cols[], const PetscScalar v[], InsertMode ins)
1066d975228cSstefano_zampini {
1067d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
1068d975228cSstefano_zampini   PetscScalar        *vals = (PetscScalar *)v;
1069d975228cSstefano_zampini   PetscScalar        sscr[MATHYPRE_SCRATCH];
10707d968826Sstefano_zampini   HYPRE_Int          cscr[2][MATHYPRE_SCRATCH];
10717d968826Sstefano_zampini   HYPRE_Int          i,nzc;
1072d975228cSstefano_zampini   PetscErrorCode     ierr;
1073d975228cSstefano_zampini 
1074d975228cSstefano_zampini   PetscFunctionBegin;
1075d975228cSstefano_zampini   for (i=0,nzc=0;i<nc;i++) {
1076d975228cSstefano_zampini     if (cols[i] >= 0) {
1077d975228cSstefano_zampini       cscr[0][nzc  ] = cols[i];
1078d975228cSstefano_zampini       cscr[1][nzc++] = i;
1079d975228cSstefano_zampini     }
1080d975228cSstefano_zampini   }
1081d975228cSstefano_zampini   if (!nzc) PetscFunctionReturn(0);
1082d975228cSstefano_zampini 
1083d975228cSstefano_zampini   if (ins == ADD_VALUES) {
1084d975228cSstefano_zampini     for (i=0;i<nr;i++) {
1085e3977e59Sstefano_zampini       if (rows[i] >= 0 && nzc) {
1086d975228cSstefano_zampini         PetscInt j;
1087d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
10887d968826Sstefano_zampini         PetscStackCallStandard(HYPRE_IJMatrixAddToValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
1089d975228cSstefano_zampini       }
1090d975228cSstefano_zampini       vals += nc;
1091d975228cSstefano_zampini     }
1092d975228cSstefano_zampini   } else { /* INSERT_VALUES */
1093d975228cSstefano_zampini #if defined(PETSC_USE_DEBUG)
1094d975228cSstefano_zampini     /* Insert values cannot be used to insert offproc entries */
1095d975228cSstefano_zampini     PetscInt rst,ren;
1096d975228cSstefano_zampini     ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
1097d975228cSstefano_zampini     for (i=0;i<nr;i++)
1098d975228cSstefano_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");
1099d975228cSstefano_zampini #endif
1100d975228cSstefano_zampini     for (i=0;i<nr;i++) {
1101e3977e59Sstefano_zampini       if (rows[i] >= 0 && nzc) {
1102d975228cSstefano_zampini         PetscInt j;
1103d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
11047d968826Sstefano_zampini         PetscStackCallStandard(HYPRE_IJMatrixSetValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
1105d975228cSstefano_zampini       }
1106d975228cSstefano_zampini       vals += nc;
1107d975228cSstefano_zampini     }
1108d975228cSstefano_zampini   }
1109d975228cSstefano_zampini   PetscFunctionReturn(0);
1110d975228cSstefano_zampini }
1111d975228cSstefano_zampini 
1112d975228cSstefano_zampini #undef __FUNCT__
1113d975228cSstefano_zampini #define __FUNCT__ "MatHYPRESetPreallocation_HYPRE"
1114d975228cSstefano_zampini static PetscErrorCode MatHYPRESetPreallocation_HYPRE(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1115d975228cSstefano_zampini {
1116d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
11177d968826Sstefano_zampini   HYPRE_Int          *hdnnz,*honnz;
111806a29025Sstefano_zampini   PetscInt           i,rs,re,cs,ce,bs;
1119d975228cSstefano_zampini   PetscMPIInt        size;
1120d975228cSstefano_zampini   PetscErrorCode     ierr;
1121d975228cSstefano_zampini 
1122d975228cSstefano_zampini   PetscFunctionBegin;
112306a29025Sstefano_zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1124d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1125d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1126d975228cSstefano_zampini   rs   = A->rmap->rstart;
1127d975228cSstefano_zampini   re   = A->rmap->rend;
1128d975228cSstefano_zampini   cs   = A->cmap->rstart;
1129d975228cSstefano_zampini   ce   = A->cmap->rend;
1130d975228cSstefano_zampini   if (!hA->ij) {
1131d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rs,re-1,cs,ce-1,&hA->ij));
1132d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1133d975228cSstefano_zampini   } else {
11347d968826Sstefano_zampini     HYPRE_Int hrs,hre,hcs,hce;
1135d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixGetLocalRange,(hA->ij,&hrs,&hre,&hcs,&hce));
1136d975228cSstefano_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);
1137d975228cSstefano_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);
1138d975228cSstefano_zampini   }
1139d975228cSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1140d975228cSstefano_zampini 
114106a29025Sstefano_zampini   if (dnz == PETSC_DEFAULT || dnz == PETSC_DECIDE) dnz = 10*bs;
114206a29025Sstefano_zampini   if (onz == PETSC_DEFAULT || onz == PETSC_DECIDE) onz = 10*bs;
114306a29025Sstefano_zampini 
1144d975228cSstefano_zampini   if (!dnnz) {
1145d975228cSstefano_zampini     ierr = PetscMalloc1(A->rmap->n,&hdnnz);CHKERRQ(ierr);
1146d975228cSstefano_zampini     for (i=0;i<A->rmap->n;i++) hdnnz[i] = dnz;
1147d975228cSstefano_zampini   } else {
11487d968826Sstefano_zampini     hdnnz = (HYPRE_Int*)dnnz;
1149d975228cSstefano_zampini   }
1150d975228cSstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1151d975228cSstefano_zampini   if (size > 1) {
1152d975228cSstefano_zampini     if (!onnz) {
1153d975228cSstefano_zampini       ierr = PetscMalloc1(A->rmap->n,&honnz);CHKERRQ(ierr);
1154d975228cSstefano_zampini       for (i=0;i<A->rmap->n;i++) honnz[i] = onz;
1155d975228cSstefano_zampini     } else {
11567d968826Sstefano_zampini       honnz = (HYPRE_Int*)onnz;
1157d975228cSstefano_zampini     }
1158d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(hA->ij,hdnnz,honnz));
1159d975228cSstefano_zampini   } else {
1160d975228cSstefano_zampini     honnz = NULL;
1161d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetRowSizes,(hA->ij,hdnnz));
1162d975228cSstefano_zampini   }
1163d975228cSstefano_zampini   if (!dnnz) {
1164d975228cSstefano_zampini     ierr = PetscFree(hdnnz);CHKERRQ(ierr);
1165d975228cSstefano_zampini   }
1166d975228cSstefano_zampini   if (!onnz && honnz) {
1167d975228cSstefano_zampini     ierr = PetscFree(honnz);CHKERRQ(ierr);
1168d975228cSstefano_zampini   }
116906a29025Sstefano_zampini   A->preallocated = PETSC_TRUE;
1170d975228cSstefano_zampini 
1171d975228cSstefano_zampini   /* SetDiagOffdSizes sets hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0 */
1172d975228cSstefano_zampini   {
1173d975228cSstefano_zampini     hypre_AuxParCSRMatrix *aux_matrix;
1174d975228cSstefano_zampini     aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij);
1175d975228cSstefano_zampini     hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1;
1176d975228cSstefano_zampini   }
1177d975228cSstefano_zampini   PetscFunctionReturn(0);
1178d975228cSstefano_zampini }
1179d975228cSstefano_zampini 
1180d975228cSstefano_zampini /*@C
1181d975228cSstefano_zampini    MatHYPRESetPreallocation - Preallocates memory for a sparse parallel matrix in HYPRE IJ format
1182d975228cSstefano_zampini 
1183d975228cSstefano_zampini    Collective on Mat
1184d975228cSstefano_zampini 
1185d975228cSstefano_zampini    Input Parameters:
1186d975228cSstefano_zampini +  A - the matrix
1187d975228cSstefano_zampini .  dnz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1188d975228cSstefano_zampini           (same value is used for all local rows)
1189d975228cSstefano_zampini .  dnnz - array containing the number of nonzeros in the various rows of the
1190d975228cSstefano_zampini           DIAGONAL portion of the local submatrix (possibly different for each row)
1191d975228cSstefano_zampini           or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure.
1192d975228cSstefano_zampini           The size of this array is equal to the number of local rows, i.e 'm'.
1193d975228cSstefano_zampini           For matrices that will be factored, you must leave room for (and set)
1194d975228cSstefano_zampini           the diagonal entry even if it is zero.
1195d975228cSstefano_zampini .  onz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1196d975228cSstefano_zampini           submatrix (same value is used for all local rows).
1197d975228cSstefano_zampini -  onnz - array containing the number of nonzeros in the various rows of the
1198d975228cSstefano_zampini           OFF-DIAGONAL portion of the local submatrix (possibly different for
1199d975228cSstefano_zampini           each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero
1200d975228cSstefano_zampini           structure. The size of this array is equal to the number
1201d975228cSstefano_zampini           of local rows, i.e 'm'.
1202d975228cSstefano_zampini 
1203d975228cSstefano_zampini    Notes: If the *nnz parameter is given then the *nz parameter is ignored; for sequential matrices, onz and onnz are ignored.
1204d975228cSstefano_zampini 
1205d975228cSstefano_zampini    Level: intermediate
1206d975228cSstefano_zampini 
1207d975228cSstefano_zampini .keywords: matrix, aij, compressed row, sparse, parallel
1208d975228cSstefano_zampini 
1209d975228cSstefano_zampini .seealso: MatCreate(), MatMPIAIJSetPreallocation, MATHYPRE
1210d975228cSstefano_zampini @*/
1211d975228cSstefano_zampini #undef __FUNCT__
1212d975228cSstefano_zampini #define __FUNCT__ "MatHYPRESetPreallocation"
1213d975228cSstefano_zampini PetscErrorCode MatHYPRESetPreallocation(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1214d975228cSstefano_zampini {
1215d975228cSstefano_zampini   PetscErrorCode ierr;
1216d975228cSstefano_zampini 
1217d975228cSstefano_zampini   PetscFunctionBegin;
1218d975228cSstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1219d975228cSstefano_zampini   PetscValidType(A,1);
1220d975228cSstefano_zampini   ierr = PetscTryMethod(A,"MatHYPRESetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(A,dnz,dnnz,onz,onnz));CHKERRQ(ierr);
1221d975228cSstefano_zampini   PetscFunctionReturn(0);
1222d975228cSstefano_zampini }
1223d975228cSstefano_zampini 
1224225daaf8SStefano Zampini /*
1225225daaf8SStefano Zampini    MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix
1226225daaf8SStefano Zampini 
1227225daaf8SStefano Zampini    Collective
1228225daaf8SStefano Zampini 
1229225daaf8SStefano Zampini    Input Parameters:
1230225daaf8SStefano Zampini +  vparcsr  - the pointer to the hypre_ParCSRMatrix
1231bb4689ddSStefano Zampini .  mtype    - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported.
1232225daaf8SStefano Zampini -  copymode - PETSc copying options
1233225daaf8SStefano Zampini 
1234225daaf8SStefano Zampini    Output Parameter:
1235225daaf8SStefano Zampini .  A  - the matrix
1236225daaf8SStefano Zampini 
1237225daaf8SStefano Zampini    Level: intermediate
1238225daaf8SStefano Zampini 
1239225daaf8SStefano Zampini .seealso: MatHYPRE, PetscCopyMode
1240225daaf8SStefano Zampini */
1241978814f1SStefano Zampini #undef __FUNCT__
1242225daaf8SStefano Zampini #define __FUNCT__ "MatCreateFromParCSR"
1243dd9c0a25Sstefano_zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(hypre_ParCSRMatrix *vparcsr, MatType mtype, PetscCopyMode copymode, Mat* A)
1244978814f1SStefano Zampini {
1245225daaf8SStefano Zampini   Mat                   T;
1246978814f1SStefano Zampini   Mat_HYPRE             *hA;
124763c07aadSStefano Zampini   hypre_ParCSRMatrix    *parcsr;
1248978814f1SStefano Zampini   MPI_Comm              comm;
1249978814f1SStefano Zampini   PetscInt              rstart,rend,cstart,cend,M,N;
1250bb4689ddSStefano Zampini   PetscBool             isseqaij,ismpiaij,isaij,ishyp,isis;
1251978814f1SStefano Zampini   PetscErrorCode        ierr;
1252978814f1SStefano Zampini 
1253978814f1SStefano Zampini   PetscFunctionBegin;
1254978814f1SStefano Zampini   parcsr = (hypre_ParCSRMatrix *)vparcsr;
1255978814f1SStefano Zampini   comm   = hypre_ParCSRMatrixComm(parcsr);
1256225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
1257225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
1258225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr);
1259225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr);
12608cfe8d00SStefano Zampini   ierr   = PetscStrcmp(mtype,MATIS,&isis);CHKERRQ(ierr);
1261225daaf8SStefano Zampini   isaij  = (PetscBool)(isseqaij || ismpiaij || isaij);
1262bb4689ddSStefano 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);
1263225daaf8SStefano Zampini   if (ishyp && copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES");
1264978814f1SStefano Zampini 
1265978814f1SStefano Zampini   /* access ParCSRMatrix */
1266978814f1SStefano Zampini   rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr);
1267978814f1SStefano Zampini   rend   = hypre_ParCSRMatrixLastRowIndex(parcsr);
1268978814f1SStefano Zampini   cstart = hypre_ParCSRMatrixFirstColDiag(parcsr);
1269978814f1SStefano Zampini   cend   = hypre_ParCSRMatrixLastColDiag(parcsr);
1270978814f1SStefano Zampini   M      = hypre_ParCSRMatrixGlobalNumRows(parcsr);
1271978814f1SStefano Zampini   N      = hypre_ParCSRMatrixGlobalNumCols(parcsr);
1272978814f1SStefano Zampini 
1273fa92c42cSstefano_zampini   /* fix for empty local rows/columns */
1274fa92c42cSstefano_zampini   if (rend < rstart) rend = rstart;
1275fa92c42cSstefano_zampini   if (cend < cstart) cend = cstart;
1276fa92c42cSstefano_zampini 
1277978814f1SStefano Zampini   /* create PETSc matrix with MatHYPRE */
1278225daaf8SStefano Zampini   ierr = MatCreate(comm,&T);CHKERRQ(ierr);
1279225daaf8SStefano Zampini   ierr = MatSetSizes(T,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr);
1280225daaf8SStefano Zampini   ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr);
1281225daaf8SStefano Zampini   hA   = (Mat_HYPRE*)(T->data);
1282978814f1SStefano Zampini 
1283978814f1SStefano Zampini   /* create HYPRE_IJMatrix */
1284978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
1285978814f1SStefano Zampini 
1286978814f1SStefano Zampini   /* set ParCSR object */
1287978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1288978814f1SStefano Zampini   hypre_IJMatrixObject(hA->ij) = parcsr;
12894ec6421dSstefano_zampini   T->preallocated = PETSC_TRUE;
1290978814f1SStefano Zampini 
1291978814f1SStefano Zampini   /* set assembled flag */
1292978814f1SStefano Zampini   hypre_IJMatrixAssembleFlag(hA->ij) = 1;
1293978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1294225daaf8SStefano Zampini   if (ishyp) {
12956d2a658fSstefano_zampini     PetscMPIInt myid = 0;
12966d2a658fSstefano_zampini 
12976d2a658fSstefano_zampini     /* make sure we always have row_starts and col_starts available */
12986d2a658fSstefano_zampini     if (HYPRE_AssumedPartitionCheck()) {
12996d2a658fSstefano_zampini       ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
13006d2a658fSstefano_zampini     }
13016d2a658fSstefano_zampini     if (!hypre_ParCSRMatrixOwnsColStarts(parcsr)) {
13026d2a658fSstefano_zampini       PetscLayout map;
13036d2a658fSstefano_zampini 
13046d2a658fSstefano_zampini       ierr = MatGetLayouts(T,NULL,&map);CHKERRQ(ierr);
13056d2a658fSstefano_zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
13061c92d2d0Sstefano_zampini       hypre_ParCSRMatrixColStarts(parcsr) = (HYPRE_Int*)(map->range + myid);
13076d2a658fSstefano_zampini     }
13086d2a658fSstefano_zampini     if (!hypre_ParCSRMatrixOwnsRowStarts(parcsr)) {
13096d2a658fSstefano_zampini       PetscLayout map;
13106d2a658fSstefano_zampini 
13116d2a658fSstefano_zampini       ierr = MatGetLayouts(T,&map,NULL);CHKERRQ(ierr);
13126d2a658fSstefano_zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
13131c92d2d0Sstefano_zampini       hypre_ParCSRMatrixRowStarts(parcsr) = (HYPRE_Int*)(map->range + myid);
13146d2a658fSstefano_zampini     }
1315978814f1SStefano Zampini     /* prevent from freeing the pointer */
1316978814f1SStefano Zampini     if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE;
1317225daaf8SStefano Zampini     *A   = T;
1318978814f1SStefano Zampini     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1319978814f1SStefano Zampini     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1320bb4689ddSStefano Zampini   } else if (isaij) {
1321bb4689ddSStefano Zampini     if (copymode != PETSC_OWN_POINTER) {
1322225daaf8SStefano Zampini       /* prevent from freeing the pointer */
1323225daaf8SStefano Zampini       hA->inner_free = PETSC_FALSE;
1324225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
1325225daaf8SStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
1326225daaf8SStefano Zampini     } else { /* AIJ return type with PETSC_OWN_POINTER */
1327225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1328225daaf8SStefano Zampini       *A   = T;
1329225daaf8SStefano Zampini     }
1330bb4689ddSStefano Zampini   } else if (isis) {
13318cfe8d00SStefano Zampini     ierr = MatConvert_HYPRE_IS(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
13328cfe8d00SStefano Zampini     if (copymode != PETSC_OWN_POINTER) hA->inner_free = PETSC_FALSE;
1333bb4689ddSStefano Zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
1334bb4689ddSStefano Zampini   }
1335978814f1SStefano Zampini   PetscFunctionReturn(0);
1336978814f1SStefano Zampini }
1337978814f1SStefano Zampini 
133863c07aadSStefano Zampini #undef __FUNCT__
1339dd9c0a25Sstefano_zampini #define __FUNCT__ "MatHYPREGetParCSR_HYPRE"
1340dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR_HYPRE(Mat A, hypre_ParCSRMatrix **parcsr)
1341dd9c0a25Sstefano_zampini {
1342dd9c0a25Sstefano_zampini   Mat_HYPRE*            hA = (Mat_HYPRE*)A->data;
1343dd9c0a25Sstefano_zampini   HYPRE_Int             type;
1344dd9c0a25Sstefano_zampini   PetscErrorCode        ierr;
1345dd9c0a25Sstefano_zampini 
1346dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1347dd9c0a25Sstefano_zampini   if (!hA->ij) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"HYPRE_IJMatrix not present");
1348dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
1349dd9c0a25Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"HYPRE_IJMatrix is not of type HYPRE_PARCSR");
1350dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)parcsr));
1351dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1352dd9c0a25Sstefano_zampini }
1353dd9c0a25Sstefano_zampini 
1354dd9c0a25Sstefano_zampini /*
1355dd9c0a25Sstefano_zampini    MatHYPREGetParCSR - Gets the pointer to the ParCSR matrix
1356dd9c0a25Sstefano_zampini 
1357dd9c0a25Sstefano_zampini    Not collective
1358dd9c0a25Sstefano_zampini 
1359dd9c0a25Sstefano_zampini    Input Parameters:
1360dd9c0a25Sstefano_zampini +  A  - the MATHYPRE object
1361dd9c0a25Sstefano_zampini 
1362dd9c0a25Sstefano_zampini    Output Parameter:
1363dd9c0a25Sstefano_zampini .  parcsr  - the pointer to the hypre_ParCSRMatrix
1364dd9c0a25Sstefano_zampini 
1365dd9c0a25Sstefano_zampini    Level: intermediate
1366dd9c0a25Sstefano_zampini 
1367dd9c0a25Sstefano_zampini .seealso: MatHYPRE, PetscCopyMode
1368dd9c0a25Sstefano_zampini */
1369dd9c0a25Sstefano_zampini #undef __FUNCT__
1370dd9c0a25Sstefano_zampini #define __FUNCT__ "MatHYPREGetParCSR"
1371dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR(Mat A, hypre_ParCSRMatrix **parcsr)
1372dd9c0a25Sstefano_zampini {
1373dd9c0a25Sstefano_zampini   PetscErrorCode ierr;
1374dd9c0a25Sstefano_zampini 
1375dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1376dd9c0a25Sstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1377dd9c0a25Sstefano_zampini   PetscValidType(A,1);
1378dd9c0a25Sstefano_zampini   ierr = PetscUseMethod(A,"MatHYPREGetParCSR_C",(Mat,hypre_ParCSRMatrix**),(A,parcsr));CHKERRQ(ierr);
1379dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1380dd9c0a25Sstefano_zampini }
1381dd9c0a25Sstefano_zampini 
1382dd9c0a25Sstefano_zampini #undef __FUNCT__
138363c07aadSStefano Zampini #define __FUNCT__ "MatCreate_HYPRE"
138463c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
138563c07aadSStefano Zampini {
138663c07aadSStefano Zampini   Mat_HYPRE      *hB;
138763c07aadSStefano Zampini   PetscErrorCode ierr;
138863c07aadSStefano Zampini 
138963c07aadSStefano Zampini   PetscFunctionBegin;
139063c07aadSStefano Zampini   ierr           = PetscNewLog(B,&hB);CHKERRQ(ierr);
1391978814f1SStefano Zampini   hB->inner_free = PETSC_TRUE;
1392978814f1SStefano Zampini 
139363c07aadSStefano Zampini   B->data       = (void*)hB;
139463c07aadSStefano Zampini   B->rmap->bs   = 1;
139563c07aadSStefano Zampini   B->assembled  = PETSC_FALSE;
139663c07aadSStefano Zampini 
139763c07aadSStefano Zampini   B->ops->mult          = MatMult_HYPRE;
139863c07aadSStefano Zampini   B->ops->multtranspose = MatMultTranspose_HYPRE;
139963c07aadSStefano Zampini   B->ops->setup         = MatSetUp_HYPRE;
140063c07aadSStefano Zampini   B->ops->destroy       = MatDestroy_HYPRE;
140163c07aadSStefano Zampini   B->ops->assemblyend   = MatAssemblyEnd_HYPRE;
1402cd8bc7baSStefano Zampini   B->ops->ptap          = MatPtAP_HYPRE_HYPRE;
1403d501dc42Sstefano_zampini   B->ops->matmult       = MatMatMult_HYPRE_HYPRE;
1404d975228cSstefano_zampini   B->ops->setvalues     = MatSetValues_HYPRE;
140563c07aadSStefano Zampini 
140663c07aadSStefano Zampini   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
140763c07aadSStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
140863c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
14092df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr);
1410c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1411c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1412d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPRESetPreallocation_C",MatHYPRESetPreallocation_HYPRE);CHKERRQ(ierr);
1413dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPREGetParCSR_C",MatHYPREGetParCSR_HYPRE);CHKERRQ(ierr);
141463c07aadSStefano Zampini   PetscFunctionReturn(0);
141563c07aadSStefano Zampini }
141663c07aadSStefano Zampini 
1417225daaf8SStefano Zampini #undef __FUNCT__
1418225daaf8SStefano Zampini #define __FUNCT__ "hypre_array_destroy"
1419225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void *ptr)
1420225daaf8SStefano Zampini {
1421225daaf8SStefano Zampini    PetscFunctionBegin;
1422225daaf8SStefano Zampini    hypre_TFree(ptr);
1423225daaf8SStefano Zampini    PetscFunctionReturn(0);
1424225daaf8SStefano Zampini }
1425225daaf8SStefano Zampini 
142663c07aadSStefano Zampini #if 0
142763c07aadSStefano Zampini /*
142863c07aadSStefano Zampini     Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format
142963c07aadSStefano Zampini 
143063c07aadSStefano Zampini     This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first
143163c07aadSStefano Zampini     which will corrupt the PETSc data structure if we did this. Need a work around to this problem.
143263c07aadSStefano Zampini */
143363c07aadSStefano Zampini #include <_hypre_IJ_mv.h>
143463c07aadSStefano Zampini #include <HYPRE_IJ_mv.h>
143563c07aadSStefano Zampini 
143663c07aadSStefano Zampini #undef __FUNCT__
143763c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixLink"
143863c07aadSStefano Zampini PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij)
143963c07aadSStefano Zampini {
144063c07aadSStefano Zampini   PetscErrorCode        ierr;
144163c07aadSStefano Zampini   PetscInt              rstart,rend,cstart,cend;
144263c07aadSStefano Zampini   PetscBool             flg;
144363c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
144463c07aadSStefano Zampini 
144563c07aadSStefano Zampini   PetscFunctionBegin;
144663c07aadSStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
144763c07aadSStefano Zampini   PetscValidType(A,1);
144863c07aadSStefano Zampini   PetscValidPointer(ij,2);
144963c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
145063c07aadSStefano Zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices");
145163c07aadSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
145263c07aadSStefano Zampini 
145363c07aadSStefano Zampini   ierr   = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
145463c07aadSStefano Zampini   rstart = A->rmap->rstart;
145563c07aadSStefano Zampini   rend   = A->rmap->rend;
145663c07aadSStefano Zampini   cstart = A->cmap->rstart;
145763c07aadSStefano Zampini   cend   = A->cmap->rend;
145863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij));
145963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
146063c07aadSStefano Zampini 
146163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij));
146263c07aadSStefano Zampini   PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij));
146363c07aadSStefano Zampini 
146463c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
146563c07aadSStefano Zampini 
146663c07aadSStefano Zampini   /* this is the Hack part where we monkey directly with the hypre datastructures */
146763c07aadSStefano Zampini 
146863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij));
146963c07aadSStefano Zampini   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
147063c07aadSStefano Zampini   PetscFunctionReturn(0);
147163c07aadSStefano Zampini }
147263c07aadSStefano Zampini #endif
1473