xref: /petsc/src/mat/impls/hypre/mhypre.c (revision bb4689dd300edcc1bdc6e337af69bbd99e9531fb)
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 
1563c07aadSStefano Zampini #include <petscsys.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;
8463c07aadSStefano Zampini   ierr   = MatSetUp(A);CHKERRQ(ierr);
8563c07aadSStefano Zampini   rstart = A->rmap->rstart;
8663c07aadSStefano Zampini   rend   = A->rmap->rend;
8763c07aadSStefano Zampini   cstart = A->cmap->rstart;
8863c07aadSStefano Zampini   cend   = A->cmap->rend;
8963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
9063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
9163c07aadSStefano Zampini   {
9263c07aadSStefano Zampini     PetscBool      same;
9363c07aadSStefano Zampini     Mat            A_d,A_o;
9463c07aadSStefano Zampini     const PetscInt *colmap;
9563c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
9663c07aadSStefano Zampini     if (same) {
9763c07aadSStefano Zampini       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
9863c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
9963c07aadSStefano Zampini       PetscFunctionReturn(0);
10063c07aadSStefano Zampini     }
10163c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
10263c07aadSStefano Zampini     if (same) {
10363c07aadSStefano Zampini       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
10463c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
10563c07aadSStefano Zampini       PetscFunctionReturn(0);
10663c07aadSStefano Zampini     }
10763c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
10863c07aadSStefano Zampini     if (same) {
10963c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
11063c07aadSStefano Zampini       PetscFunctionReturn(0);
11163c07aadSStefano Zampini     }
11263c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
11363c07aadSStefano Zampini     if (same) {
11463c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
11563c07aadSStefano Zampini       PetscFunctionReturn(0);
11663c07aadSStefano Zampini     }
11763c07aadSStefano Zampini   }
11863c07aadSStefano Zampini   PetscFunctionReturn(0);
11963c07aadSStefano Zampini }
12063c07aadSStefano Zampini 
12163c07aadSStefano Zampini #undef __FUNCT__
12263c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixCopy"
12363c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij)
12463c07aadSStefano Zampini {
12563c07aadSStefano Zampini   PetscErrorCode    ierr;
12663c07aadSStefano Zampini   PetscInt          i,rstart,rend,ncols,nr,nc;
12763c07aadSStefano Zampini   const PetscScalar *values;
12863c07aadSStefano Zampini   const PetscInt    *cols;
12963c07aadSStefano Zampini   PetscBool         flg;
13063c07aadSStefano Zampini 
13163c07aadSStefano Zampini   PetscFunctionBegin;
13263c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
13363c07aadSStefano Zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
13463c07aadSStefano Zampini   if (flg && nr == nc) {
13563c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
13663c07aadSStefano Zampini     PetscFunctionReturn(0);
13763c07aadSStefano Zampini   }
13863c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
13963c07aadSStefano Zampini   if (flg) {
14063c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
14163c07aadSStefano Zampini     PetscFunctionReturn(0);
14263c07aadSStefano Zampini   }
14363c07aadSStefano Zampini 
14463c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
14563c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
14663c07aadSStefano Zampini   for (i=rstart; i<rend; i++) {
14763c07aadSStefano Zampini     ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
14863c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
14963c07aadSStefano Zampini     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
15063c07aadSStefano Zampini   }
15163c07aadSStefano Zampini   PetscFunctionReturn(0);
15263c07aadSStefano Zampini }
15363c07aadSStefano Zampini 
15463c07aadSStefano Zampini #undef __FUNCT__
15563c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ"
15663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij)
15763c07aadSStefano Zampini {
15863c07aadSStefano Zampini   PetscErrorCode        ierr;
15963c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag = (Mat_SeqAIJ*)A->data;
16058968eb6SStefano Zampini   HYPRE_Int             type;
16163c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
16263c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
16363c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag;
16463c07aadSStefano Zampini 
16563c07aadSStefano Zampini   PetscFunctionBegin;
16663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
167ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
168ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
169ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
17063c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
17163c07aadSStefano Zampini   /*
17263c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
17363c07aadSStefano Zampini   */
17463c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));
17563c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));
17663c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
177ea9daf28SStefano Zampini 
178ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
17963c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
18063c07aadSStefano Zampini   PetscFunctionReturn(0);
18163c07aadSStefano Zampini }
18263c07aadSStefano Zampini 
18363c07aadSStefano Zampini #undef __FUNCT__
18463c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ"
18563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
18663c07aadSStefano Zampini {
18763c07aadSStefano Zampini   PetscErrorCode        ierr;
18863c07aadSStefano Zampini   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
18963c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag,*poffd;
19063c07aadSStefano Zampini   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
19158968eb6SStefano Zampini   HYPRE_Int             type;
19263c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
19363c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
19463c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag,*hoffd;
19563c07aadSStefano Zampini 
19663c07aadSStefano Zampini   PetscFunctionBegin;
19763c07aadSStefano Zampini   pdiag = (Mat_SeqAIJ*) pA->A->data;
19863c07aadSStefano Zampini   poffd = (Mat_SeqAIJ*) pA->B->data;
19963c07aadSStefano Zampini   /* cstart is only valid for square MPIAIJ layed out in the usual way */
20063c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
20163c07aadSStefano Zampini 
20263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
203ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
204ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
205ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
20663c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
20763c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
20863c07aadSStefano Zampini 
20963c07aadSStefano Zampini   /*
21063c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
21163c07aadSStefano Zampini   */
21263c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
21363c07aadSStefano Zampini   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
21463c07aadSStefano Zampini   jj  = (PetscInt*)hdiag->j;
21563c07aadSStefano Zampini   pjj = (PetscInt*)pdiag->j;
21663c07aadSStefano Zampini   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
21763c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
21863c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
21963c07aadSStefano Zampini   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
22063c07aadSStefano Zampini      If we hacked a hypre a bit more we might be able to avoid this step */
22163c07aadSStefano Zampini   jj  = (PetscInt*) hoffd->j;
22263c07aadSStefano Zampini   pjj = (PetscInt*) poffd->j;
22363c07aadSStefano Zampini   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
22463c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
22563c07aadSStefano Zampini 
226ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
22763c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
22863c07aadSStefano Zampini   PetscFunctionReturn(0);
22963c07aadSStefano Zampini }
23063c07aadSStefano Zampini 
23163c07aadSStefano Zampini #undef __FUNCT__
2322df22349SStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_IS"
2332df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B)
2342df22349SStefano Zampini {
2352df22349SStefano Zampini   Mat_HYPRE*             mhA = (Mat_HYPRE*)(A->data);
2362df22349SStefano Zampini   Mat                    lA;
2372df22349SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
2382df22349SStefano Zampini   IS                     is;
2392df22349SStefano Zampini   hypre_ParCSRMatrix     *hA;
2402df22349SStefano Zampini   hypre_CSRMatrix        *hdiag,*hoffd;
2412df22349SStefano Zampini   MPI_Comm               comm;
2422df22349SStefano Zampini   PetscScalar            *hdd,*hod,*aa,*data;
2432df22349SStefano Zampini   PetscInt               *col_map_offd,*hdi,*hdj,*hoi,*hoj;
2442df22349SStefano Zampini   PetscInt               *ii,*jj,*iptr,*jptr;
2452df22349SStefano Zampini   PetscInt               cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N;
24658968eb6SStefano Zampini   HYPRE_Int              type;
2472df22349SStefano Zampini   PetscErrorCode         ierr;
2482df22349SStefano Zampini 
2492df22349SStefano Zampini   PetscFunctionBegin;
250a1787963SStefano Zampini   comm = PetscObjectComm((PetscObject)A);
2512df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type));
2522df22349SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
2532df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA));
2542df22349SStefano Zampini   M     = hypre_ParCSRMatrixGlobalNumRows(hA);
2552df22349SStefano Zampini   N     = hypre_ParCSRMatrixGlobalNumCols(hA);
2562df22349SStefano Zampini   str   = hypre_ParCSRMatrixFirstRowIndex(hA);
2572df22349SStefano Zampini   stc   = hypre_ParCSRMatrixFirstColDiag(hA);
2582df22349SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(hA);
2592df22349SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(hA);
2602df22349SStefano Zampini   dr    = hypre_CSRMatrixNumRows(hdiag);
2612df22349SStefano Zampini   dc    = hypre_CSRMatrixNumCols(hdiag);
2622df22349SStefano Zampini   nnz   = hypre_CSRMatrixNumNonzeros(hdiag);
2632df22349SStefano Zampini   hdi   = hypre_CSRMatrixI(hdiag);
2642df22349SStefano Zampini   hdj   = hypre_CSRMatrixJ(hdiag);
2652df22349SStefano Zampini   hdd   = hypre_CSRMatrixData(hdiag);
2662df22349SStefano Zampini   oc    = hypre_CSRMatrixNumCols(hoffd);
2672df22349SStefano Zampini   nnz  += hypre_CSRMatrixNumNonzeros(hoffd);
2682df22349SStefano Zampini   hoi   = hypre_CSRMatrixI(hoffd);
2692df22349SStefano Zampini   hoj   = hypre_CSRMatrixJ(hoffd);
2702df22349SStefano Zampini   hod   = hypre_CSRMatrixData(hoffd);
2712df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
2722df22349SStefano Zampini     PetscInt *aux;
2732df22349SStefano Zampini 
2742df22349SStefano Zampini     /* generate l2g maps for rows and cols */
2752df22349SStefano Zampini     ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr);
2762df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
2772df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2782df22349SStefano Zampini     col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
2792df22349SStefano Zampini     ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr);
2802df22349SStefano Zampini     for (i=0; i<dc; i++) aux[i] = i+stc;
2812df22349SStefano Zampini     for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i];
2822df22349SStefano Zampini     ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
2832df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
2842df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2852df22349SStefano Zampini     /* create MATIS object */
2862df22349SStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
2872df22349SStefano Zampini     ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr);
2882df22349SStefano Zampini     ierr = MatSetType(*B,MATIS);CHKERRQ(ierr);
2892df22349SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr);
2902df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
2912df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
2922df22349SStefano Zampini 
2932df22349SStefano Zampini     /* allocate CSR for local matrix */
2942df22349SStefano Zampini     ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr);
2952df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr);
2962df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr);
2972df22349SStefano Zampini   } else {
2982df22349SStefano Zampini     PetscInt  nr;
2992df22349SStefano Zampini     PetscBool done;
3002df22349SStefano Zampini     ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr);
3012df22349SStefano Zampini     ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr);
3022df22349SStefano 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);
3032df22349SStefano 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);
3042df22349SStefano Zampini     ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr);
3052df22349SStefano Zampini   }
3062df22349SStefano Zampini   /* merge local matrices */
3072df22349SStefano Zampini   ii   = iptr;
3082df22349SStefano Zampini   jj   = jptr;
3092df22349SStefano Zampini   aa   = data;
3102df22349SStefano Zampini   *ii  = *(hdi++) + *(hoi++);
3112df22349SStefano Zampini   for (jd=0,jo=0,cum=0; *ii<nnz; cum++) {
3122df22349SStefano Zampini     PetscScalar *aold = aa;
3132df22349SStefano Zampini     PetscInt    *jold = jj,nc = jd+jo;
3142df22349SStefano Zampini     for (; jd<*hdi; jd++) { *jj++ = *hdj++;      *aa++ = *hdd++; }
3152df22349SStefano Zampini     for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
3162df22349SStefano Zampini     *(++ii) = *(hdi++) + *(hoi++);
3172df22349SStefano Zampini     ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr);
3182df22349SStefano Zampini   }
3192df22349SStefano Zampini   for (; cum<dr; cum++) *(++ii) = nnz;
3202df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
321a033916dSStefano Zampini     Mat_SeqAIJ* a;
322a033916dSStefano Zampini 
3232df22349SStefano Zampini     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr);
3242df22349SStefano Zampini     ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr);
325a033916dSStefano Zampini     /* hack SeqAIJ */
326a033916dSStefano Zampini     a          = (Mat_SeqAIJ*)(lA->data);
327a033916dSStefano Zampini     a->free_a  = PETSC_TRUE;
328a033916dSStefano Zampini     a->free_ij = PETSC_TRUE;
3292df22349SStefano Zampini     ierr = MatDestroy(&lA);CHKERRQ(ierr);
3302df22349SStefano Zampini   }
3312df22349SStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3322df22349SStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3332df22349SStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
3342df22349SStefano Zampini     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
3352df22349SStefano Zampini   }
3362df22349SStefano Zampini   PetscFunctionReturn(0);
3372df22349SStefano Zampini }
3382df22349SStefano Zampini 
3392df22349SStefano Zampini #undef __FUNCT__
34063c07aadSStefano Zampini #define __FUNCT__ "MatConvert_AIJ_HYPRE"
34163c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
34263c07aadSStefano Zampini {
34363c07aadSStefano Zampini   Mat_HYPRE      *hB;
34463c07aadSStefano Zampini   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
34563c07aadSStefano Zampini   PetscErrorCode ierr;
34663c07aadSStefano Zampini 
34763c07aadSStefano Zampini   PetscFunctionBegin;
34863c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
34963c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
35063c07aadSStefano Zampini     /* always destroy the old matrix and create a new memory;
35163c07aadSStefano Zampini        hope this does not churn the memory too much. The problem
35263c07aadSStefano Zampini        is I do not know if it is possible to put the matrix back to
35363c07aadSStefano Zampini        its initial state so that we can directly copy the values
35463c07aadSStefano Zampini        the second time through. */
35563c07aadSStefano Zampini     hB = (Mat_HYPRE*)((*B)->data);
35663c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
35763c07aadSStefano Zampini   } else {
35863c07aadSStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
35963c07aadSStefano Zampini     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
36063c07aadSStefano Zampini     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
36163c07aadSStefano Zampini     ierr = MatSetUp(*B);CHKERRQ(ierr);
36263c07aadSStefano Zampini     hB   = (Mat_HYPRE*)((*B)->data);
36363c07aadSStefano Zampini   }
36463c07aadSStefano Zampini   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
36563c07aadSStefano Zampini   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
36663c07aadSStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
36763c07aadSStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
36863c07aadSStefano Zampini   PetscFunctionReturn(0);
36963c07aadSStefano Zampini }
37063c07aadSStefano Zampini 
37163c07aadSStefano Zampini #undef __FUNCT__
37263c07aadSStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_AIJ"
373ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
37463c07aadSStefano Zampini {
37563c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
37663c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
37763c07aadSStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
37863c07aadSStefano Zampini   MPI_Comm           comm;
37963c07aadSStefano Zampini   PetscScalar        *da,*oa,*aptr;
38063c07aadSStefano Zampini   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
38163c07aadSStefano Zampini   PetscInt           i,dnnz,onnz,m,n;
38258968eb6SStefano Zampini   HYPRE_Int          type;
38363c07aadSStefano Zampini   PetscMPIInt        size;
38463c07aadSStefano Zampini   PetscErrorCode     ierr;
38563c07aadSStefano Zampini 
38663c07aadSStefano Zampini   PetscFunctionBegin;
38763c07aadSStefano Zampini   comm = PetscObjectComm((PetscObject)A);
38863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
38963c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
39063c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
39163c07aadSStefano Zampini     PetscBool ismpiaij,isseqaij;
39263c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
39363c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
39463c07aadSStefano Zampini     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
39563c07aadSStefano Zampini   }
39663c07aadSStefano Zampini   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
39763c07aadSStefano Zampini 
39863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
39963c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(parcsr);
40063c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(parcsr);
40163c07aadSStefano Zampini   m     = hypre_CSRMatrixNumRows(hdiag);
40263c07aadSStefano Zampini   n     = hypre_CSRMatrixNumCols(hdiag);
40363c07aadSStefano Zampini   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
40463c07aadSStefano Zampini   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
405225daaf8SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
40663c07aadSStefano Zampini     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
40763c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
40863c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
409225daaf8SStefano Zampini   } else if (reuse == MAT_REUSE_MATRIX) {
41063c07aadSStefano Zampini     PetscInt  nr;
41163c07aadSStefano Zampini     PetscBool done;
41263c07aadSStefano Zampini     if (size > 1) {
41363c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
41463c07aadSStefano Zampini 
41563c07aadSStefano Zampini       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
41663c07aadSStefano 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);
41763c07aadSStefano 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);
41863c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
41963c07aadSStefano Zampini     } else {
42063c07aadSStefano Zampini       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
42163c07aadSStefano Zampini       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
42263c07aadSStefano 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);
42363c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
42463c07aadSStefano Zampini     }
425225daaf8SStefano Zampini   } else { /* MAT_INPLACE_MATRIX */
426225daaf8SStefano Zampini     dii = hypre_CSRMatrixI(hdiag);
427225daaf8SStefano Zampini     djj = hypre_CSRMatrixJ(hdiag);
428225daaf8SStefano Zampini     da  = hypre_CSRMatrixData(hdiag);
42963c07aadSStefano Zampini   }
43063c07aadSStefano Zampini   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
43163c07aadSStefano Zampini   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
43263c07aadSStefano Zampini   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
43363c07aadSStefano Zampini   iptr = djj;
43463c07aadSStefano Zampini   aptr = da;
43563c07aadSStefano Zampini   for (i=0; i<m; i++) {
43663c07aadSStefano Zampini     PetscInt nc = dii[i+1]-dii[i];
43763c07aadSStefano Zampini     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
43863c07aadSStefano Zampini     iptr += nc;
43963c07aadSStefano Zampini     aptr += nc;
44063c07aadSStefano Zampini   }
44163c07aadSStefano Zampini   if (size > 1) {
44263c07aadSStefano Zampini     PetscInt *offdj,*coffd;
44363c07aadSStefano Zampini 
444225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
44563c07aadSStefano Zampini       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
44663c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
44763c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
448225daaf8SStefano Zampini     } else if (reuse == MAT_REUSE_MATRIX) {
44963c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
45063c07aadSStefano Zampini       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
45163c07aadSStefano Zampini       PetscBool  done;
45263c07aadSStefano Zampini 
45363c07aadSStefano Zampini       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
45463c07aadSStefano 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);
45563c07aadSStefano 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);
45663c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
457225daaf8SStefano Zampini     } else { /* MAT_INPLACE_MATRIX */
458225daaf8SStefano Zampini       oii = hypre_CSRMatrixI(hoffd);
459225daaf8SStefano Zampini       ojj = hypre_CSRMatrixJ(hoffd);
460225daaf8SStefano Zampini       oa  = hypre_CSRMatrixData(hoffd);
46163c07aadSStefano Zampini     }
46263c07aadSStefano Zampini     ierr  = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
46363c07aadSStefano Zampini     offdj = hypre_CSRMatrixJ(hoffd);
46463c07aadSStefano Zampini     coffd = hypre_ParCSRMatrixColMapOffd(parcsr);
46563c07aadSStefano Zampini     for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]];
46663c07aadSStefano Zampini     ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr);
46763c07aadSStefano Zampini     iptr = ojj;
46863c07aadSStefano Zampini     aptr = oa;
46963c07aadSStefano Zampini     for (i=0; i<m; i++) {
47063c07aadSStefano Zampini        PetscInt nc = oii[i+1]-oii[i];
47163c07aadSStefano Zampini        ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
47263c07aadSStefano Zampini        iptr += nc;
47363c07aadSStefano Zampini        aptr += nc;
47463c07aadSStefano Zampini     }
475225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
47663c07aadSStefano Zampini       Mat_MPIAIJ *b;
47763c07aadSStefano Zampini       Mat_SeqAIJ *d,*o;
478225daaf8SStefano Zampini 
47963c07aadSStefano Zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,
48063c07aadSStefano Zampini                                             djj,da,oii,ojj,oa,B);CHKERRQ(ierr);
48163c07aadSStefano Zampini       /* hack MPIAIJ */
48263c07aadSStefano Zampini       b          = (Mat_MPIAIJ*)((*B)->data);
48363c07aadSStefano Zampini       d          = (Mat_SeqAIJ*)b->A->data;
48463c07aadSStefano Zampini       o          = (Mat_SeqAIJ*)b->B->data;
48563c07aadSStefano Zampini       d->free_a  = PETSC_TRUE;
48663c07aadSStefano Zampini       d->free_ij = PETSC_TRUE;
48763c07aadSStefano Zampini       o->free_a  = PETSC_TRUE;
48863c07aadSStefano Zampini       o->free_ij = PETSC_TRUE;
489225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
490225daaf8SStefano Zampini       Mat T;
491225daaf8SStefano Zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,
492225daaf8SStefano Zampini                                             djj,da,oii,ojj,oa,&T);CHKERRQ(ierr);
493225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
494225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
495225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
496225daaf8SStefano Zampini       hypre_CSRMatrixI(hoffd)    = NULL;
497225daaf8SStefano Zampini       hypre_CSRMatrixJ(hoffd)    = NULL;
498225daaf8SStefano Zampini       hypre_CSRMatrixData(hoffd) = NULL;
499225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
50063c07aadSStefano Zampini     }
501225daaf8SStefano Zampini   } else {
502225daaf8SStefano Zampini     oii  = NULL;
503225daaf8SStefano Zampini     ojj  = NULL;
504225daaf8SStefano Zampini     oa   = NULL;
505225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
50663c07aadSStefano Zampini       Mat_SeqAIJ* b;
50763c07aadSStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
50863c07aadSStefano Zampini       /* hack SeqAIJ */
50963c07aadSStefano Zampini       b          = (Mat_SeqAIJ*)((*B)->data);
51063c07aadSStefano Zampini       b->free_a  = PETSC_TRUE;
51163c07aadSStefano Zampini       b->free_ij = PETSC_TRUE;
512225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
513225daaf8SStefano Zampini       Mat T;
514225daaf8SStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,&T);CHKERRQ(ierr);
515225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
516225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
517225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
518225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
51963c07aadSStefano Zampini     }
520225daaf8SStefano Zampini   }
521225daaf8SStefano Zampini 
522225daaf8SStefano Zampini   /* we have to use hypre_Tfree to free the arrays */
52363c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
524225daaf8SStefano Zampini     void *ptrs[6] = {dii,djj,da,oii,ojj,oa};
525225daaf8SStefano Zampini     const char *names[6] = {"_hypre_csr_dii",
526225daaf8SStefano Zampini                             "_hypre_csr_djj",
527225daaf8SStefano Zampini                             "_hypre_csr_da",
528225daaf8SStefano Zampini                             "_hypre_csr_oii",
529225daaf8SStefano Zampini                             "_hypre_csr_ojj",
530225daaf8SStefano Zampini                             "_hypre_csr_oa"};
531225daaf8SStefano Zampini     for (i=0; i<6; i++) {
532225daaf8SStefano Zampini       PetscContainer c;
533225daaf8SStefano Zampini 
534225daaf8SStefano Zampini       ierr = PetscContainerCreate(comm,&c);CHKERRQ(ierr);
535225daaf8SStefano Zampini       ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr);
536225daaf8SStefano Zampini       ierr = PetscContainerSetUserDestroy(c,hypre_array_destroy);CHKERRQ(ierr);
537225daaf8SStefano Zampini       ierr = PetscObjectCompose((PetscObject)(*B),names[i],(PetscObject)c);CHKERRQ(ierr);
538225daaf8SStefano Zampini       ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
539225daaf8SStefano Zampini     }
54063c07aadSStefano Zampini   }
54163c07aadSStefano Zampini   PetscFunctionReturn(0);
54263c07aadSStefano Zampini }
54363c07aadSStefano Zampini 
54463c07aadSStefano Zampini #undef __FUNCT__
545c1a070e6SStefano Zampini #define __FUNCT__ "MatPtAP_AIJ_HYPRE"
546c1a070e6SStefano Zampini static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
547c1a070e6SStefano Zampini {
548c1a070e6SStefano Zampini   Mat_HYPRE          *hP = (Mat_HYPRE*)P->data;
549c1a070e6SStefano Zampini   hypre_ParCSRMatrix *tA,*Pparcsr,*ptapparcsr;
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   HYPRE_Int          type,P_owns_col_starts;
554c1a070e6SStefano Zampini   PetscBool          ismpiaij,isseqaij;
555c1a070e6SStefano Zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
556225daaf8SStefano Zampini   char               mtype[256];
557c1a070e6SStefano Zampini   PetscErrorCode     ierr;
558c1a070e6SStefano Zampini 
559c1a070e6SStefano Zampini   PetscFunctionBegin;
560c1a070e6SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
561c1a070e6SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
562c1a070e6SStefano Zampini   if (scall == MAT_REUSE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_REUSE_MATRIX");
563c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
564c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
565c1a070e6SStefano Zampini   if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type);
566225daaf8SStefano Zampini 
567225daaf8SStefano Zampini   /* It looks like we don't need to have the diagonal entries
568225daaf8SStefano Zampini      ordered first in the rows of the diagonal part
569225daaf8SStefano Zampini      for boomerAMGBuildCoarseOperator to work */
570c1a070e6SStefano Zampini   if (ismpiaij) {
571c1a070e6SStefano Zampini     Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data);
572c1a070e6SStefano Zampini 
573c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)a->A->data;
574c1a070e6SStefano Zampini     offd   = (Mat_SeqAIJ*)a->B->data;
575c1a070e6SStefano Zampini     garray = a->garray;
576c1a070e6SStefano Zampini     noffd  = a->B->cmap->N;
577c1a070e6SStefano Zampini     dnnz   = diag->nz;
578c1a070e6SStefano Zampini     onnz   = offd->nz;
579c1a070e6SStefano Zampini   } else {
580c1a070e6SStefano Zampini     diag    = (Mat_SeqAIJ*)A->data;
581c1a070e6SStefano Zampini     offd    = NULL;
582c1a070e6SStefano Zampini     garray  = NULL;
583c1a070e6SStefano Zampini     noffd   = 0;
584c1a070e6SStefano Zampini     dnnz    = diag->nz;
585c1a070e6SStefano Zampini     onnz    = 0;
586c1a070e6SStefano Zampini   }
587225daaf8SStefano Zampini 
588c1a070e6SStefano Zampini   /* create a temporary ParCSR */
589c1a070e6SStefano Zampini   if (HYPRE_AssumedPartitionCheck()) {
590c1a070e6SStefano Zampini    PetscMPIInt myid;
591c1a070e6SStefano Zampini 
592c1a070e6SStefano Zampini    ierr       = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
593c1a070e6SStefano Zampini    row_starts = A->rmap->range + myid;
594c1a070e6SStefano Zampini    col_starts = A->cmap->range + myid;
595c1a070e6SStefano Zampini   } else {
596c1a070e6SStefano Zampini    row_starts = A->rmap->range;
597c1a070e6SStefano Zampini    col_starts = A->cmap->range;
598c1a070e6SStefano Zampini   }
599c1a070e6SStefano Zampini   tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,row_starts,col_starts,noffd,dnnz,onnz);
600c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(tA,0);
601c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(tA,0);
602c1a070e6SStefano Zampini 
603225daaf8SStefano Zampini   /* set diagonal part */
604c1a070e6SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(tA);
605c1a070e6SStefano Zampini   hypre_CSRMatrixI(hdiag)           = diag->i;
606c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hdiag)           = diag->j;
607c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = diag->a;
608c1a070e6SStefano Zampini   hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz;
609c1a070e6SStefano Zampini   hypre_CSRMatrixSetRownnz(hdiag);
610c1a070e6SStefano Zampini   hypre_CSRMatrixSetDataOwner(hdiag,0);
611c1a070e6SStefano Zampini 
612225daaf8SStefano Zampini   /* set offdiagonal part */
613c1a070e6SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(tA);
614c1a070e6SStefano Zampini   if (offd) {
615c1a070e6SStefano Zampini     hypre_CSRMatrixI(hoffd)           = offd->i;
616c1a070e6SStefano Zampini     hypre_CSRMatrixJ(hoffd)           = offd->j;
617c1a070e6SStefano Zampini     hypre_CSRMatrixData(hoffd)        = offd->a;
618c1a070e6SStefano Zampini     hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz;
619c1a070e6SStefano Zampini     hypre_CSRMatrixSetRownnz(hoffd);
620c1a070e6SStefano Zampini     hypre_CSRMatrixSetDataOwner(hoffd,0);
621c1a070e6SStefano Zampini     hypre_ParCSRMatrixSetNumNonzeros(tA);
622c1a070e6SStefano Zampini     hypre_ParCSRMatrixColMapOffd(tA) = garray;
623c1a070e6SStefano Zampini   }
624c1a070e6SStefano Zampini 
625c1a070e6SStefano Zampini   /* call RAP from BoomerAMG */
626c1a070e6SStefano Zampini   /* Warning: hypre_BoomerAMGBuildCoarseOperator steals the col_starts
627c1a070e6SStefano Zampini      from Pparcsr (even if it does not own them)! */
628c1a070e6SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
629c1a070e6SStefano Zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(Pparcsr);
630c1a070e6SStefano Zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
631c1a070e6SStefano Zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(Pparcsr,tA,Pparcsr,&ptapparcsr));
632c1a070e6SStefano Zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
633c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(ptapparcsr,0);
634c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(ptapparcsr,0);
635c1a070e6SStefano Zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(Pparcsr, 1);
636c1a070e6SStefano Zampini 
637c1a070e6SStefano Zampini   /* set pointers to NULL before destroying tA */
638c1a070e6SStefano Zampini   hypre_CSRMatrixI(hdiag)          = NULL;
639c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hdiag)          = NULL;
640c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)       = NULL;
641c1a070e6SStefano Zampini   hypre_CSRMatrixI(hoffd)          = NULL;
642c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hoffd)          = NULL;
643c1a070e6SStefano Zampini   hypre_CSRMatrixData(hoffd)       = NULL;
644c1a070e6SStefano Zampini   hypre_ParCSRMatrixColMapOffd(tA) = NULL;
645c1a070e6SStefano Zampini   hypre_ParCSRMatrixDestroy(tA);
646225daaf8SStefano Zampini 
647225daaf8SStefano Zampini   /* create C depending on mtype */
648225daaf8SStefano Zampini   sprintf(mtype,MATAIJ);
649225daaf8SStefano Zampini   ierr = PetscOptionsGetString(((PetscObject)A)->options,((PetscObject)A)->prefix,"-matptap_hypre_outtype",mtype,256,NULL);CHKERRQ(ierr);
650225daaf8SStefano Zampini   ierr = MatCreateFromParCSR(ptapparcsr,mtype,PETSC_OWN_POINTER,C);CHKERRQ(ierr);
651c1a070e6SStefano Zampini   PetscFunctionReturn(0);
652c1a070e6SStefano Zampini }
653c1a070e6SStefano Zampini 
654c1a070e6SStefano Zampini #undef __FUNCT__
655cd8bc7baSStefano Zampini #define __FUNCT__ "MatPtAP_HYPRE_HYPRE"
656cd8bc7baSStefano Zampini static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
657cd8bc7baSStefano Zampini {
658cd8bc7baSStefano Zampini   hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr;
659cd8bc7baSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data, *hP = (Mat_HYPRE*)P->data;
660c1a070e6SStefano Zampini   HYPRE_Int          type,P_owns_col_starts;
661cd8bc7baSStefano Zampini   PetscErrorCode     ierr;
662cd8bc7baSStefano Zampini 
663cd8bc7baSStefano Zampini   PetscFunctionBegin;
664cd8bc7baSStefano Zampini   if (scall == MAT_REUSE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Unsupported MAT_REUSE_MATRIX");
665cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
666cd8bc7baSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
667cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
668cd8bc7baSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
669cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
670cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
671c1a070e6SStefano Zampini 
672c1a070e6SStefano Zampini   /* call RAP from BoomerAMG */
673c1a070e6SStefano Zampini   /* Warning: hypre_BoomerAMGBuildCoarseOperator steals the col_starts
674c1a070e6SStefano Zampini      from Pparcsr (even if it does not own them)! */
675c1a070e6SStefano Zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(Pparcsr);
676c1a070e6SStefano Zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
677cd8bc7baSStefano Zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr));
678cd8bc7baSStefano Zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
679c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(ptapparcsr,0);
680c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(ptapparcsr,0);
681c1a070e6SStefano Zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(Pparcsr, 1);
682c1a070e6SStefano Zampini 
683c1a070e6SStefano Zampini   /* create MatHYPRE */
684225daaf8SStefano Zampini   ierr = MatCreateFromParCSR(ptapparcsr,MATHYPRE,PETSC_OWN_POINTER,C);CHKERRQ(ierr);
685cd8bc7baSStefano Zampini   PetscFunctionReturn(0);
686cd8bc7baSStefano Zampini }
687cd8bc7baSStefano Zampini 
688cd8bc7baSStefano Zampini #undef __FUNCT__
68963c07aadSStefano Zampini #define __FUNCT__ "MatMultTranspose_HYPRE"
690ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
69163c07aadSStefano Zampini {
69263c07aadSStefano Zampini   PetscErrorCode ierr;
69363c07aadSStefano Zampini 
69463c07aadSStefano Zampini   PetscFunctionBegin;
69563c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
69663c07aadSStefano Zampini   PetscFunctionReturn(0);
69763c07aadSStefano Zampini }
69863c07aadSStefano Zampini 
69963c07aadSStefano Zampini #undef __FUNCT__
70063c07aadSStefano Zampini #define __FUNCT__ "MatMult_HYPRE"
701ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
70263c07aadSStefano Zampini {
70363c07aadSStefano Zampini   PetscErrorCode ierr;
70463c07aadSStefano Zampini 
70563c07aadSStefano Zampini   PetscFunctionBegin;
70663c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
70763c07aadSStefano Zampini   PetscFunctionReturn(0);
70863c07aadSStefano Zampini }
70963c07aadSStefano Zampini 
71063c07aadSStefano Zampini #undef __FUNCT__
71163c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_MultKernel_Private"
712ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
71363c07aadSStefano Zampini {
71463c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
71563c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
71663c07aadSStefano Zampini   hypre_ParVector    *hx,*hy;
71763c07aadSStefano Zampini   PetscScalar        *ax,*ay,*sax,*say;
71863c07aadSStefano Zampini   PetscErrorCode     ierr;
71963c07aadSStefano Zampini 
72063c07aadSStefano Zampini   PetscFunctionBegin;
72163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
72263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
72363c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
72463c07aadSStefano Zampini   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
72563c07aadSStefano Zampini   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
72663c07aadSStefano Zampini   if (trans) {
72758968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ay,say);
72858968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax);
72963c07aadSStefano Zampini     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
73058968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,say,ay);
73158968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax);
73263c07aadSStefano Zampini   } else {
73358968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax);
73458968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ay,say);
73563c07aadSStefano Zampini     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
73658968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax);
73758968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,say,ay);
73863c07aadSStefano Zampini   }
73963c07aadSStefano Zampini   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
74063c07aadSStefano Zampini   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
74163c07aadSStefano Zampini   PetscFunctionReturn(0);
74263c07aadSStefano Zampini }
74363c07aadSStefano Zampini 
74463c07aadSStefano Zampini #undef __FUNCT__
74563c07aadSStefano Zampini #define __FUNCT__ "MatDestroy_HYPRE"
746ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A)
74763c07aadSStefano Zampini {
74863c07aadSStefano Zampini   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
74963c07aadSStefano Zampini   PetscErrorCode ierr;
75063c07aadSStefano Zampini 
75163c07aadSStefano Zampini   PetscFunctionBegin;
75263c07aadSStefano Zampini   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
75363c07aadSStefano Zampini   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
754978814f1SStefano Zampini   if (hA->ij) {
755978814f1SStefano Zampini     if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL;
756978814f1SStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
757978814f1SStefano Zampini   }
75863c07aadSStefano Zampini   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); }
75963c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
7602df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr);
761c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr);
762c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr);
76363c07aadSStefano Zampini   ierr = PetscFree(A->data);CHKERRQ(ierr);
76463c07aadSStefano Zampini   PetscFunctionReturn(0);
76563c07aadSStefano Zampini }
76663c07aadSStefano Zampini 
76763c07aadSStefano Zampini #undef __FUNCT__
76863c07aadSStefano Zampini #define __FUNCT__ "MatSetUp_HYPRE"
769ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A)
77063c07aadSStefano Zampini {
77163c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
77263c07aadSStefano Zampini   Vec                x,b;
77363c07aadSStefano Zampini   PetscErrorCode     ierr;
77463c07aadSStefano Zampini 
77563c07aadSStefano Zampini   PetscFunctionBegin;
77663c07aadSStefano Zampini   if (hA->x) PetscFunctionReturn(0);
77763c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
77863c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
77963c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
78063c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
78163c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
78263c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
78363c07aadSStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
78463c07aadSStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
78563c07aadSStefano Zampini   PetscFunctionReturn(0);
78663c07aadSStefano Zampini }
78763c07aadSStefano Zampini 
78863c07aadSStefano Zampini #undef __FUNCT__
78963c07aadSStefano Zampini #define __FUNCT__ "MatAssemblyEnd_HYPRE"
790ea9daf28SStefano Zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
79163c07aadSStefano Zampini {
79263c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
79363c07aadSStefano Zampini   PetscErrorCode     ierr;
79463c07aadSStefano Zampini 
79563c07aadSStefano Zampini   PetscFunctionBegin;
79663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
79763c07aadSStefano Zampini   PetscFunctionReturn(0);
79863c07aadSStefano Zampini }
79963c07aadSStefano Zampini 
800225daaf8SStefano Zampini /*
801225daaf8SStefano Zampini    MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix
802225daaf8SStefano Zampini 
803225daaf8SStefano Zampini    Collective
804225daaf8SStefano Zampini 
805225daaf8SStefano Zampini    Input Parameters:
806225daaf8SStefano Zampini +  vparcsr  - the pointer to the hypre_ParCSRMatrix
807*bb4689ddSStefano Zampini .  mtype    - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported.
808225daaf8SStefano Zampini -  copymode - PETSc copying options
809225daaf8SStefano Zampini 
810225daaf8SStefano Zampini    Output Parameter:
811225daaf8SStefano Zampini .  A  - the matrix
812225daaf8SStefano Zampini 
813225daaf8SStefano Zampini    Level: intermediate
814225daaf8SStefano Zampini 
815*bb4689ddSStefano Zampini    Notes: Not all the combinations between copymode and mtype are currently supported.
816225daaf8SStefano Zampini 
817225daaf8SStefano Zampini .seealso: MatHYPRE, PetscCopyMode
818225daaf8SStefano Zampini */
819978814f1SStefano Zampini #undef __FUNCT__
820225daaf8SStefano Zampini #define __FUNCT__ "MatCreateFromParCSR"
821225daaf8SStefano Zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(void *vparcsr, MatType mtype, PetscCopyMode copymode, Mat* A)
822978814f1SStefano Zampini {
823225daaf8SStefano Zampini   Mat                   T;
824978814f1SStefano Zampini   Mat_HYPRE             *hA;
82563c07aadSStefano Zampini   hypre_ParCSRMatrix    *parcsr;
826978814f1SStefano Zampini   MPI_Comm              comm;
827978814f1SStefano Zampini   PetscInt              rstart,rend,cstart,cend,M,N;
828*bb4689ddSStefano Zampini   PetscBool             isseqaij,ismpiaij,isaij,ishyp,isis;
829978814f1SStefano Zampini   PetscErrorCode        ierr;
830978814f1SStefano Zampini 
831978814f1SStefano Zampini   PetscFunctionBegin;
832978814f1SStefano Zampini   parcsr = (hypre_ParCSRMatrix *)vparcsr;
833978814f1SStefano Zampini   comm   = hypre_ParCSRMatrixComm(parcsr);
834225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
835225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
836225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr);
837225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr);
838*bb4689ddSStefano Zampini   ierr   = PetscStrcmp(mtype,MATHYPRE,&isis);CHKERRQ(ierr);
839225daaf8SStefano Zampini   isaij  = (PetscBool)(isseqaij || ismpiaij || isaij);
840*bb4689ddSStefano 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);
841225daaf8SStefano Zampini   if (ishyp && copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES");
842*bb4689ddSStefano Zampini   if (isis  && copymode == PETSC_OWN_POINTER) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_OWN_POINTER");
843978814f1SStefano Zampini 
844978814f1SStefano Zampini   /* access ParCSRMatrix */
845978814f1SStefano Zampini   rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr);
846978814f1SStefano Zampini   rend   = hypre_ParCSRMatrixLastRowIndex(parcsr);
847978814f1SStefano Zampini   cstart = hypre_ParCSRMatrixFirstColDiag(parcsr);
848978814f1SStefano Zampini   cend   = hypre_ParCSRMatrixLastColDiag(parcsr);
849978814f1SStefano Zampini   M      = hypre_ParCSRMatrixGlobalNumRows(parcsr);
850978814f1SStefano Zampini   N      = hypre_ParCSRMatrixGlobalNumCols(parcsr);
851978814f1SStefano Zampini 
852978814f1SStefano Zampini   /* create PETSc matrix with MatHYPRE */
853225daaf8SStefano Zampini   ierr = MatCreate(comm,&T);CHKERRQ(ierr);
854225daaf8SStefano Zampini   ierr = MatSetSizes(T,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr);
855225daaf8SStefano Zampini   ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr);
856225daaf8SStefano Zampini   ierr = MatSetUp(T);CHKERRQ(ierr);
857225daaf8SStefano Zampini   hA   = (Mat_HYPRE*)(T->data);
858978814f1SStefano Zampini 
859978814f1SStefano Zampini   /* create HYPRE_IJMatrix */
860978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
861978814f1SStefano Zampini 
862978814f1SStefano Zampini   /* set ParCSR object */
863978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
864978814f1SStefano Zampini   hypre_IJMatrixObject(hA->ij) = parcsr;
865978814f1SStefano Zampini 
866978814f1SStefano Zampini   /* set assembled flag */
867978814f1SStefano Zampini   hypre_IJMatrixAssembleFlag(hA->ij) = 1;
868978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
869225daaf8SStefano Zampini   if (ishyp) {
870978814f1SStefano Zampini     /* prevent from freeing the pointer */
871978814f1SStefano Zampini     if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE;
872225daaf8SStefano Zampini     *A = T;
873978814f1SStefano Zampini     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
874978814f1SStefano Zampini     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
875*bb4689ddSStefano Zampini   } else if (isaij) {
876*bb4689ddSStefano Zampini     if (copymode != PETSC_OWN_POINTER) {
877225daaf8SStefano Zampini       /* prevent from freeing the pointer */
878225daaf8SStefano Zampini       hA->inner_free = PETSC_FALSE;
879225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
880225daaf8SStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
881225daaf8SStefano Zampini     } else { /* AIJ return type with PETSC_OWN_POINTER */
882225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
883225daaf8SStefano Zampini       *A   = T;
884225daaf8SStefano Zampini     }
885*bb4689ddSStefano Zampini   } else if (isis) {
886*bb4689ddSStefano Zampini     ierr = MatConvert_HYPRE_AIJ(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
887*bb4689ddSStefano Zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
888*bb4689ddSStefano Zampini   }
889978814f1SStefano Zampini   PetscFunctionReturn(0);
890978814f1SStefano Zampini }
891978814f1SStefano Zampini 
89263c07aadSStefano Zampini #undef __FUNCT__
89363c07aadSStefano Zampini #define __FUNCT__ "MatCreate_HYPRE"
89463c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
89563c07aadSStefano Zampini {
89663c07aadSStefano Zampini   Mat_HYPRE      *hB;
89763c07aadSStefano Zampini   PetscErrorCode ierr;
89863c07aadSStefano Zampini 
89963c07aadSStefano Zampini   PetscFunctionBegin;
90063c07aadSStefano Zampini   ierr           = PetscNewLog(B,&hB);CHKERRQ(ierr);
901978814f1SStefano Zampini   hB->inner_free = PETSC_TRUE;
902978814f1SStefano Zampini 
90363c07aadSStefano Zampini   B->data       = (void*)hB;
90463c07aadSStefano Zampini   B->rmap->bs   = 1;
90563c07aadSStefano Zampini   B->assembled  = PETSC_FALSE;
90663c07aadSStefano Zampini 
90763c07aadSStefano Zampini   B->ops->mult          = MatMult_HYPRE;
90863c07aadSStefano Zampini   B->ops->multtranspose = MatMultTranspose_HYPRE;
90963c07aadSStefano Zampini   B->ops->setup         = MatSetUp_HYPRE;
91063c07aadSStefano Zampini   B->ops->destroy       = MatDestroy_HYPRE;
91163c07aadSStefano Zampini   B->ops->assemblyend   = MatAssemblyEnd_HYPRE;
912cd8bc7baSStefano Zampini   B->ops->ptap          = MatPtAP_HYPRE_HYPRE;
91363c07aadSStefano Zampini 
91463c07aadSStefano Zampini   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
91563c07aadSStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
91663c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
9172df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr);
918c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
919c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
92063c07aadSStefano Zampini   PetscFunctionReturn(0);
92163c07aadSStefano Zampini }
92263c07aadSStefano Zampini 
923225daaf8SStefano Zampini #undef __FUNCT__
924225daaf8SStefano Zampini #define __FUNCT__ "hypre_array_destroy"
925225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void *ptr)
926225daaf8SStefano Zampini {
927225daaf8SStefano Zampini    PetscFunctionBegin;
928225daaf8SStefano Zampini    hypre_TFree(ptr);
929225daaf8SStefano Zampini    PetscFunctionReturn(0);
930225daaf8SStefano Zampini }
931225daaf8SStefano Zampini 
93263c07aadSStefano Zampini #if 0
93363c07aadSStefano Zampini /*
93463c07aadSStefano Zampini     Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format
93563c07aadSStefano Zampini 
93663c07aadSStefano Zampini     This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first
93763c07aadSStefano Zampini     which will corrupt the PETSc data structure if we did this. Need a work around to this problem.
93863c07aadSStefano Zampini */
93963c07aadSStefano Zampini #include <_hypre_IJ_mv.h>
94063c07aadSStefano Zampini #include <HYPRE_IJ_mv.h>
94163c07aadSStefano Zampini 
94263c07aadSStefano Zampini #undef __FUNCT__
94363c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixLink"
94463c07aadSStefano Zampini PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij)
94563c07aadSStefano Zampini {
94663c07aadSStefano Zampini   PetscErrorCode        ierr;
94763c07aadSStefano Zampini   PetscInt              rstart,rend,cstart,cend;
94863c07aadSStefano Zampini   PetscBool             flg;
94963c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
95063c07aadSStefano Zampini 
95163c07aadSStefano Zampini   PetscFunctionBegin;
95263c07aadSStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
95363c07aadSStefano Zampini   PetscValidType(A,1);
95463c07aadSStefano Zampini   PetscValidPointer(ij,2);
95563c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
95663c07aadSStefano Zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices");
95763c07aadSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
95863c07aadSStefano Zampini 
95963c07aadSStefano Zampini   ierr   = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
96063c07aadSStefano Zampini   rstart = A->rmap->rstart;
96163c07aadSStefano Zampini   rend   = A->rmap->rend;
96263c07aadSStefano Zampini   cstart = A->cmap->rstart;
96363c07aadSStefano Zampini   cend   = A->cmap->rend;
96463c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij));
96563c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
96663c07aadSStefano Zampini 
96763c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij));
96863c07aadSStefano Zampini   PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij));
96963c07aadSStefano Zampini 
97063c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
97163c07aadSStefano Zampini 
97263c07aadSStefano Zampini   /* this is the Hack part where we monkey directly with the hypre datastructures */
97363c07aadSStefano Zampini 
97463c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij));
97563c07aadSStefano Zampini   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
97663c07aadSStefano Zampini   PetscFunctionReturn(0);
97763c07aadSStefano Zampini }
97863c07aadSStefano Zampini #endif
979