xref: /petsc/src/mat/impls/hypre/mhypre.c (revision e3977e59a21d22c294746612581e7b3969ff4c18)
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);
149*e3977e59Sstefano_zampini     if (ncols) {
15063c07aadSStefano Zampini       PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
151*e3977e59Sstefano_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:
642613e5ff0Sstefano_zampini    the resulting ParCSR will not own the column and row starts */
643613e5ff0Sstefano_zampini #undef __FUNCT__
644613e5ff0Sstefano_zampini #define __FUNCT__ "MatHYPRE_ParCSR_RAP"
645613e5ff0Sstefano_zampini static PetscErrorCode MatHYPRE_ParCSR_RAP(hypre_ParCSRMatrix *hR, hypre_ParCSRMatrix *hA,
646613e5ff0Sstefano_zampini                                           hypre_ParCSRMatrix *hP, hypre_ParCSRMatrix **hRAP)
647613e5ff0Sstefano_zampini {
648613e5ff0Sstefano_zampini   PetscErrorCode ierr;
649613e5ff0Sstefano_zampini   HYPRE_Int      P_owns_col_starts,R_owns_row_starts;
650613e5ff0Sstefano_zampini 
651613e5ff0Sstefano_zampini   PetscFunctionBegin;
652613e5ff0Sstefano_zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(hP);
653613e5ff0Sstefano_zampini   R_owns_row_starts = hypre_ParCSRMatrixOwnsRowStarts(hR);
654613e5ff0Sstefano_zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(hR,hA,hP,hRAP));
655613e5ff0Sstefano_zampini   PetscStackCallStandard(hypre_ParCSRMatrixSetNumNonzeros,(*hRAP));
656613e5ff0Sstefano_zampini   /* hypre_BoomerAMGBuildCoarseOperator steals the col_starts from P and the row_starts from R */
657613e5ff0Sstefano_zampini   hypre_ParCSRMatrixSetRowStartsOwner(*hRAP,0);
658613e5ff0Sstefano_zampini   hypre_ParCSRMatrixSetColStartsOwner(*hRAP,0);
659613e5ff0Sstefano_zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(hP,1);
660613e5ff0Sstefano_zampini   if (R_owns_row_starts) hypre_ParCSRMatrixSetRowStartsOwner(hR,1);
661613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
662613e5ff0Sstefano_zampini }
663613e5ff0Sstefano_zampini 
664613e5ff0Sstefano_zampini #undef __FUNCT__
6656f231fbdSstefano_zampini #define __FUNCT__ "MatPtAPNumeric_AIJ_AIJ_wHYPRE"
6666f231fbdSstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat P,Mat C)
667613e5ff0Sstefano_zampini {
6686f231fbdSstefano_zampini   Mat                B;
669613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *hA,*hP,*hPtAP;
670613e5ff0Sstefano_zampini   PetscErrorCode     ierr;
671613e5ff0Sstefano_zampini 
672613e5ff0Sstefano_zampini   PetscFunctionBegin;
673613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
674613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(P,&hP);CHKERRQ(ierr);
675613e5ff0Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(hP,hA,hP,&hPtAP);CHKERRQ(ierr);
6766f231fbdSstefano_zampini   ierr = MatCreateFromParCSR(hPtAP,MATAIJ,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
6776f231fbdSstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
678613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
679613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(P,&hP);CHKERRQ(ierr);
6806f231fbdSstefano_zampini   PetscFunctionReturn(0);
6816f231fbdSstefano_zampini }
6826f231fbdSstefano_zampini 
6836f231fbdSstefano_zampini #undef __FUNCT__
6846f231fbdSstefano_zampini #define __FUNCT__ "MatPtAPSymbolic_AIJ_AIJ_wHYPRE"
6856f231fbdSstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat P,PetscReal fill,Mat *C)
6866f231fbdSstefano_zampini {
6876f231fbdSstefano_zampini   PetscErrorCode     ierr;
6886f231fbdSstefano_zampini   PetscFunctionBegin;
6896f231fbdSstefano_zampini   ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
6906f231fbdSstefano_zampini   ierr                   = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
6916f231fbdSstefano_zampini   (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_AIJ_wHYPRE;
692613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
693613e5ff0Sstefano_zampini }
694613e5ff0Sstefano_zampini 
695613e5ff0Sstefano_zampini #undef __FUNCT__
6964cc28894Sstefano_zampini #define __FUNCT__ "MatPtAPNumeric_AIJ_HYPRE"
6974cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_HYPRE(Mat A,Mat P,Mat C)
698613e5ff0Sstefano_zampini {
6994cc28894Sstefano_zampini   Mat                B;
7004cc28894Sstefano_zampini   Mat_HYPRE          *hP;
701613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *hA = NULL,*Pparcsr,*ptapparcsr;
702613e5ff0Sstefano_zampini   HYPRE_Int          type;
703613e5ff0Sstefano_zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
7044cc28894Sstefano_zampini   PetscBool          ishypre;
705613e5ff0Sstefano_zampini   PetscErrorCode     ierr;
706613e5ff0Sstefano_zampini 
707613e5ff0Sstefano_zampini   PetscFunctionBegin;
7084cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
7094cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(comm,PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
7104cc28894Sstefano_zampini   hP = (Mat_HYPRE*)P->data;
711613e5ff0Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
712613e5ff0Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
713613e5ff0Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
714613e5ff0Sstefano_zampini 
715613e5ff0Sstefano_zampini   /* It looks like we don't need to have the diagonal entries
716613e5ff0Sstefano_zampini      ordered first in the rows of the diagonal part
717613e5ff0Sstefano_zampini      for boomerAMGBuildCoarseOperator to work */
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 
809cd8bc7baSStefano Zampini #undef __FUNCT__
8105e5acdf2Sstefano_zampini #define __FUNCT__ "MatMatMultNumeric_AIJ_AIJ_wHYPRE"
8115e5acdf2Sstefano_zampini static PetscErrorCode MatMatMultNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C)
8125e5acdf2Sstefano_zampini {
8135e5acdf2Sstefano_zampini   Mat                D;
8145e5acdf2Sstefano_zampini   hypre_ParCSRMatrix *hA,*hB,*hAB;
8155e5acdf2Sstefano_zampini   PetscErrorCode     ierr;
8165e5acdf2Sstefano_zampini 
8175e5acdf2Sstefano_zampini   PetscFunctionBegin;
8185e5acdf2Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
8195e5acdf2Sstefano_zampini   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
8205e5acdf2Sstefano_zampini   PetscStackPush("hypre_ParMatmul");
8215e5acdf2Sstefano_zampini   hAB  = hypre_ParMatmul(hA,hB);
8225e5acdf2Sstefano_zampini   PetscStackPop;
8235e5acdf2Sstefano_zampini   ierr = MatCreateFromParCSR(hAB,MATAIJ,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
8245e5acdf2Sstefano_zampini   ierr = MatHeaderMerge(C,&D);CHKERRQ(ierr);
8255e5acdf2Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
8265e5acdf2Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
8275e5acdf2Sstefano_zampini   PetscFunctionReturn(0);
8285e5acdf2Sstefano_zampini }
8295e5acdf2Sstefano_zampini 
8305e5acdf2Sstefano_zampini #undef __FUNCT__
8315e5acdf2Sstefano_zampini #define __FUNCT__ "MatMatMultSymbolic_AIJ_AIJ_wHYPRE"
8325e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat B,PetscReal fill,Mat *C)
8335e5acdf2Sstefano_zampini {
8345e5acdf2Sstefano_zampini   PetscErrorCode ierr;
83520e1dc0dSstefano_zampini 
8365e5acdf2Sstefano_zampini   PetscFunctionBegin;
8375e5acdf2Sstefano_zampini   ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
8385e5acdf2Sstefano_zampini   ierr                      = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
8395e5acdf2Sstefano_zampini   (*C)->ops->matmultnumeric = MatMatMultNumeric_AIJ_AIJ_wHYPRE;
8405e5acdf2Sstefano_zampini   PetscFunctionReturn(0);
8415e5acdf2Sstefano_zampini }
8425e5acdf2Sstefano_zampini 
8435e5acdf2Sstefano_zampini #undef __FUNCT__
84420e1dc0dSstefano_zampini #define __FUNCT__ "MatMatMatMultNumeric_AIJ_AIJ_AIJ_wHYPRE"
84520e1dc0dSstefano_zampini static PetscErrorCode MatMatMatMultNumeric_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,Mat D)
84620e1dc0dSstefano_zampini {
84720e1dc0dSstefano_zampini   Mat                E;
84820e1dc0dSstefano_zampini   hypre_ParCSRMatrix *hA,*hB,*hC,*hABC;
84920e1dc0dSstefano_zampini   PetscErrorCode     ierr;
85020e1dc0dSstefano_zampini 
85120e1dc0dSstefano_zampini   PetscFunctionBegin;
85220e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
85320e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
85420e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(C,&hC);CHKERRQ(ierr);
85520e1dc0dSstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(hA,hB,hC,&hABC);CHKERRQ(ierr);
85620e1dc0dSstefano_zampini   ierr = MatCreateFromParCSR(hABC,MATAIJ,PETSC_OWN_POINTER,&E);CHKERRQ(ierr);
85720e1dc0dSstefano_zampini   ierr = MatHeaderMerge(D,&E);CHKERRQ(ierr);
85820e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
85920e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
86020e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(C,&hC);CHKERRQ(ierr);
86120e1dc0dSstefano_zampini   PetscFunctionReturn(0);
86220e1dc0dSstefano_zampini }
86320e1dc0dSstefano_zampini 
86420e1dc0dSstefano_zampini #undef __FUNCT__
86520e1dc0dSstefano_zampini #define __FUNCT__ "MatMatMatMultSymbolic_AIJ_AIJ_AIJ_wHYPRE"
86620e1dc0dSstefano_zampini PETSC_INTERN PetscErrorCode MatMatMatMultSymbolic_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,PetscReal fill,Mat *D)
86720e1dc0dSstefano_zampini {
86820e1dc0dSstefano_zampini   PetscErrorCode ierr;
86920e1dc0dSstefano_zampini 
87020e1dc0dSstefano_zampini   PetscFunctionBegin;
87120e1dc0dSstefano_zampini   ierr                         = MatCreate(PetscObjectComm((PetscObject)A),D);CHKERRQ(ierr);
87220e1dc0dSstefano_zampini   ierr                         = MatSetType(*D,MATAIJ);CHKERRQ(ierr);
87320e1dc0dSstefano_zampini   (*D)->ops->matmatmultnumeric = MatMatMatMultNumeric_AIJ_AIJ_AIJ_wHYPRE;
87420e1dc0dSstefano_zampini   PetscFunctionReturn(0);
87520e1dc0dSstefano_zampini }
87620e1dc0dSstefano_zampini 
87720e1dc0dSstefano_zampini #undef __FUNCT__
87863c07aadSStefano Zampini #define __FUNCT__ "MatMultTranspose_HYPRE"
879ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
88063c07aadSStefano Zampini {
88163c07aadSStefano Zampini   PetscErrorCode ierr;
88263c07aadSStefano Zampini 
88363c07aadSStefano Zampini   PetscFunctionBegin;
88463c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
88563c07aadSStefano Zampini   PetscFunctionReturn(0);
88663c07aadSStefano Zampini }
88763c07aadSStefano Zampini 
88863c07aadSStefano Zampini #undef __FUNCT__
88963c07aadSStefano Zampini #define __FUNCT__ "MatMult_HYPRE"
890ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
89163c07aadSStefano Zampini {
89263c07aadSStefano Zampini   PetscErrorCode ierr;
89363c07aadSStefano Zampini 
89463c07aadSStefano Zampini   PetscFunctionBegin;
89563c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
89663c07aadSStefano Zampini   PetscFunctionReturn(0);
89763c07aadSStefano Zampini }
89863c07aadSStefano Zampini 
89963c07aadSStefano Zampini #undef __FUNCT__
90063c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_MultKernel_Private"
901ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
90263c07aadSStefano Zampini {
90363c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
90463c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
90563c07aadSStefano Zampini   hypre_ParVector    *hx,*hy;
90663c07aadSStefano Zampini   PetscScalar        *ax,*ay,*sax,*say;
90763c07aadSStefano Zampini   PetscErrorCode     ierr;
90863c07aadSStefano Zampini 
90963c07aadSStefano Zampini   PetscFunctionBegin;
91063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
91163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
91263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
91363c07aadSStefano Zampini   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
91463c07aadSStefano Zampini   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
91563c07aadSStefano Zampini   if (trans) {
91658968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ay,say);
91758968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax);
91863c07aadSStefano Zampini     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
91958968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,say,ay);
92058968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax);
92163c07aadSStefano Zampini   } else {
92258968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax);
92358968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ay,say);
92463c07aadSStefano Zampini     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
92558968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax);
92658968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,say,ay);
92763c07aadSStefano Zampini   }
92863c07aadSStefano Zampini   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
92963c07aadSStefano Zampini   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
93063c07aadSStefano Zampini   PetscFunctionReturn(0);
93163c07aadSStefano Zampini }
93263c07aadSStefano Zampini 
93363c07aadSStefano Zampini #undef __FUNCT__
93463c07aadSStefano Zampini #define __FUNCT__ "MatDestroy_HYPRE"
935ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A)
93663c07aadSStefano Zampini {
93763c07aadSStefano Zampini   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
93863c07aadSStefano Zampini   PetscErrorCode ierr;
93963c07aadSStefano Zampini 
94063c07aadSStefano Zampini   PetscFunctionBegin;
94163c07aadSStefano Zampini   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
94263c07aadSStefano Zampini   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
943978814f1SStefano Zampini   if (hA->ij) {
944978814f1SStefano Zampini     if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL;
945978814f1SStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
946978814f1SStefano Zampini   }
94763c07aadSStefano Zampini   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); }
94863c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
9492df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr);
950c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr);
951c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr);
952d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPRESetPreallocation_C",NULL);CHKERRQ(ierr);
953dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPREGetParCSR_C",NULL);CHKERRQ(ierr);
95463c07aadSStefano Zampini   ierr = PetscFree(A->data);CHKERRQ(ierr);
95563c07aadSStefano Zampini   PetscFunctionReturn(0);
95663c07aadSStefano Zampini }
95763c07aadSStefano Zampini 
95863c07aadSStefano Zampini #undef __FUNCT__
95963c07aadSStefano Zampini #define __FUNCT__ "MatSetUp_HYPRE"
960ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A)
96163c07aadSStefano Zampini {
9624ec6421dSstefano_zampini   PetscErrorCode ierr;
9634ec6421dSstefano_zampini 
9644ec6421dSstefano_zampini   PetscFunctionBegin;
9654ec6421dSstefano_zampini   ierr = MatHYPRESetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr);
9664ec6421dSstefano_zampini   PetscFunctionReturn(0);
9674ec6421dSstefano_zampini }
9684ec6421dSstefano_zampini 
9694ec6421dSstefano_zampini #undef __FUNCT__
9704ec6421dSstefano_zampini #define __FUNCT__ "MatAssemblyEnd_HYPRE"
9714ec6421dSstefano_zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
9724ec6421dSstefano_zampini {
97363c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
97463c07aadSStefano Zampini   Vec                x,b;
97563c07aadSStefano Zampini   PetscErrorCode     ierr;
97663c07aadSStefano Zampini 
97763c07aadSStefano Zampini   PetscFunctionBegin;
9784ec6421dSstefano_zampini   if (mode == MAT_FLUSH_ASSEMBLY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_FLUSH_ASSEMBLY currently not supported with MATHYPRE");
9794ec6421dSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
98063c07aadSStefano Zampini   if (hA->x) PetscFunctionReturn(0);
98163c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
98263c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
98363c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
98463c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
98563c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
98663c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
98763c07aadSStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
98863c07aadSStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
98963c07aadSStefano Zampini   PetscFunctionReturn(0);
99063c07aadSStefano Zampini }
99163c07aadSStefano Zampini 
992d975228cSstefano_zampini #define MATHYPRE_SCRATCH 2048
993d975228cSstefano_zampini 
994d975228cSstefano_zampini #undef __FUNCT__
995d975228cSstefano_zampini #define __FUNCT__ "MatSetValues_HYPRE"
996d975228cSstefano_zampini PetscErrorCode MatSetValues_HYPRE(Mat A, PetscInt nr, const PetscInt rows[], PetscInt nc, const PetscInt cols[], const PetscScalar v[], InsertMode ins)
997d975228cSstefano_zampini {
998d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
999d975228cSstefano_zampini   PetscScalar        *vals = (PetscScalar *)v;
1000d975228cSstefano_zampini   PetscScalar        sscr[MATHYPRE_SCRATCH];
10017d968826Sstefano_zampini   HYPRE_Int          cscr[2][MATHYPRE_SCRATCH];
10027d968826Sstefano_zampini   HYPRE_Int          i,nzc;
1003d975228cSstefano_zampini   PetscErrorCode     ierr;
1004d975228cSstefano_zampini 
1005d975228cSstefano_zampini   PetscFunctionBegin;
1006d975228cSstefano_zampini   for (i=0,nzc=0;i<nc;i++) {
1007d975228cSstefano_zampini     if (cols[i] >= 0) {
1008d975228cSstefano_zampini       cscr[0][nzc  ] = cols[i];
1009d975228cSstefano_zampini       cscr[1][nzc++] = i;
1010d975228cSstefano_zampini     }
1011d975228cSstefano_zampini   }
1012d975228cSstefano_zampini   if (!nzc) PetscFunctionReturn(0);
1013d975228cSstefano_zampini 
1014d975228cSstefano_zampini   if (ins == ADD_VALUES) {
1015d975228cSstefano_zampini     for (i=0;i<nr;i++) {
1016*e3977e59Sstefano_zampini       if (rows[i] >= 0 && nzc) {
1017d975228cSstefano_zampini         PetscInt j;
1018d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
10197d968826Sstefano_zampini         PetscStackCallStandard(HYPRE_IJMatrixAddToValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
1020d975228cSstefano_zampini       }
1021d975228cSstefano_zampini       vals += nc;
1022d975228cSstefano_zampini     }
1023d975228cSstefano_zampini   } else { /* INSERT_VALUES */
1024d975228cSstefano_zampini #if defined(PETSC_USE_DEBUG)
1025d975228cSstefano_zampini     /* Insert values cannot be used to insert offproc entries */
1026d975228cSstefano_zampini     PetscInt rst,ren;
1027d975228cSstefano_zampini     ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
1028d975228cSstefano_zampini     for (i=0;i<nr;i++)
1029d975228cSstefano_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");
1030d975228cSstefano_zampini #endif
1031d975228cSstefano_zampini     for (i=0;i<nr;i++) {
1032*e3977e59Sstefano_zampini       if (rows[i] >= 0 && nzc) {
1033d975228cSstefano_zampini         PetscInt j;
1034d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
10357d968826Sstefano_zampini         PetscStackCallStandard(HYPRE_IJMatrixSetValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
1036d975228cSstefano_zampini       }
1037d975228cSstefano_zampini       vals += nc;
1038d975228cSstefano_zampini     }
1039d975228cSstefano_zampini   }
1040d975228cSstefano_zampini   PetscFunctionReturn(0);
1041d975228cSstefano_zampini }
1042d975228cSstefano_zampini 
1043d975228cSstefano_zampini #undef __FUNCT__
1044d975228cSstefano_zampini #define __FUNCT__ "MatHYPRESetPreallocation_HYPRE"
1045d975228cSstefano_zampini static PetscErrorCode MatHYPRESetPreallocation_HYPRE(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1046d975228cSstefano_zampini {
1047d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
10487d968826Sstefano_zampini   HYPRE_Int          *hdnnz,*honnz;
104906a29025Sstefano_zampini   PetscInt           i,rs,re,cs,ce,bs;
1050d975228cSstefano_zampini   PetscMPIInt        size;
1051d975228cSstefano_zampini   PetscErrorCode     ierr;
1052d975228cSstefano_zampini 
1053d975228cSstefano_zampini   PetscFunctionBegin;
105406a29025Sstefano_zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1055d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1056d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1057d975228cSstefano_zampini   rs   = A->rmap->rstart;
1058d975228cSstefano_zampini   re   = A->rmap->rend;
1059d975228cSstefano_zampini   cs   = A->cmap->rstart;
1060d975228cSstefano_zampini   ce   = A->cmap->rend;
1061d975228cSstefano_zampini   if (!hA->ij) {
1062d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rs,re-1,cs,ce-1,&hA->ij));
1063d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1064d975228cSstefano_zampini   } else {
10657d968826Sstefano_zampini     HYPRE_Int hrs,hre,hcs,hce;
1066d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixGetLocalRange,(hA->ij,&hrs,&hre,&hcs,&hce));
1067d975228cSstefano_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);
1068d975228cSstefano_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);
1069d975228cSstefano_zampini   }
1070d975228cSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1071d975228cSstefano_zampini 
107206a29025Sstefano_zampini   if (dnz == PETSC_DEFAULT || dnz == PETSC_DECIDE) dnz = 10*bs;
107306a29025Sstefano_zampini   if (onz == PETSC_DEFAULT || onz == PETSC_DECIDE) onz = 10*bs;
107406a29025Sstefano_zampini 
1075d975228cSstefano_zampini   if (!dnnz) {
1076d975228cSstefano_zampini     ierr = PetscMalloc1(A->rmap->n,&hdnnz);CHKERRQ(ierr);
1077d975228cSstefano_zampini     for (i=0;i<A->rmap->n;i++) hdnnz[i] = dnz;
1078d975228cSstefano_zampini   } else {
10797d968826Sstefano_zampini     hdnnz = (HYPRE_Int*)dnnz;
1080d975228cSstefano_zampini   }
1081d975228cSstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1082d975228cSstefano_zampini   if (size > 1) {
1083d975228cSstefano_zampini     if (!onnz) {
1084d975228cSstefano_zampini       ierr = PetscMalloc1(A->rmap->n,&honnz);CHKERRQ(ierr);
1085d975228cSstefano_zampini       for (i=0;i<A->rmap->n;i++) honnz[i] = onz;
1086d975228cSstefano_zampini     } else {
10877d968826Sstefano_zampini       honnz = (HYPRE_Int*)onnz;
1088d975228cSstefano_zampini     }
1089d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(hA->ij,hdnnz,honnz));
1090d975228cSstefano_zampini   } else {
1091d975228cSstefano_zampini     honnz = NULL;
1092d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetRowSizes,(hA->ij,hdnnz));
1093d975228cSstefano_zampini   }
1094d975228cSstefano_zampini   if (!dnnz) {
1095d975228cSstefano_zampini     ierr = PetscFree(hdnnz);CHKERRQ(ierr);
1096d975228cSstefano_zampini   }
1097d975228cSstefano_zampini   if (!onnz && honnz) {
1098d975228cSstefano_zampini     ierr = PetscFree(honnz);CHKERRQ(ierr);
1099d975228cSstefano_zampini   }
110006a29025Sstefano_zampini   A->preallocated = PETSC_TRUE;
1101d975228cSstefano_zampini 
1102d975228cSstefano_zampini   /* SetDiagOffdSizes sets hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0 */
1103d975228cSstefano_zampini   {
1104d975228cSstefano_zampini     hypre_AuxParCSRMatrix *aux_matrix;
1105d975228cSstefano_zampini     aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij);
1106d975228cSstefano_zampini     hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1;
1107d975228cSstefano_zampini   }
1108d975228cSstefano_zampini   PetscFunctionReturn(0);
1109d975228cSstefano_zampini }
1110d975228cSstefano_zampini 
1111d975228cSstefano_zampini /*@C
1112d975228cSstefano_zampini    MatHYPRESetPreallocation - Preallocates memory for a sparse parallel matrix in HYPRE IJ format
1113d975228cSstefano_zampini 
1114d975228cSstefano_zampini    Collective on Mat
1115d975228cSstefano_zampini 
1116d975228cSstefano_zampini    Input Parameters:
1117d975228cSstefano_zampini +  A - the matrix
1118d975228cSstefano_zampini .  dnz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1119d975228cSstefano_zampini           (same value is used for all local rows)
1120d975228cSstefano_zampini .  dnnz - array containing the number of nonzeros in the various rows of the
1121d975228cSstefano_zampini           DIAGONAL portion of the local submatrix (possibly different for each row)
1122d975228cSstefano_zampini           or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure.
1123d975228cSstefano_zampini           The size of this array is equal to the number of local rows, i.e 'm'.
1124d975228cSstefano_zampini           For matrices that will be factored, you must leave room for (and set)
1125d975228cSstefano_zampini           the diagonal entry even if it is zero.
1126d975228cSstefano_zampini .  onz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1127d975228cSstefano_zampini           submatrix (same value is used for all local rows).
1128d975228cSstefano_zampini -  onnz - array containing the number of nonzeros in the various rows of the
1129d975228cSstefano_zampini           OFF-DIAGONAL portion of the local submatrix (possibly different for
1130d975228cSstefano_zampini           each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero
1131d975228cSstefano_zampini           structure. The size of this array is equal to the number
1132d975228cSstefano_zampini           of local rows, i.e 'm'.
1133d975228cSstefano_zampini 
1134d975228cSstefano_zampini    Notes: If the *nnz parameter is given then the *nz parameter is ignored; for sequential matrices, onz and onnz are ignored.
1135d975228cSstefano_zampini 
1136d975228cSstefano_zampini    Level: intermediate
1137d975228cSstefano_zampini 
1138d975228cSstefano_zampini .keywords: matrix, aij, compressed row, sparse, parallel
1139d975228cSstefano_zampini 
1140d975228cSstefano_zampini .seealso: MatCreate(), MatMPIAIJSetPreallocation, MATHYPRE
1141d975228cSstefano_zampini @*/
1142d975228cSstefano_zampini #undef __FUNCT__
1143d975228cSstefano_zampini #define __FUNCT__ "MatHYPRESetPreallocation"
1144d975228cSstefano_zampini PetscErrorCode MatHYPRESetPreallocation(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1145d975228cSstefano_zampini {
1146d975228cSstefano_zampini   PetscErrorCode ierr;
1147d975228cSstefano_zampini 
1148d975228cSstefano_zampini   PetscFunctionBegin;
1149d975228cSstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1150d975228cSstefano_zampini   PetscValidType(A,1);
1151d975228cSstefano_zampini   ierr = PetscTryMethod(A,"MatHYPRESetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(A,dnz,dnnz,onz,onnz));CHKERRQ(ierr);
1152d975228cSstefano_zampini   PetscFunctionReturn(0);
1153d975228cSstefano_zampini }
1154d975228cSstefano_zampini 
1155225daaf8SStefano Zampini /*
1156225daaf8SStefano Zampini    MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix
1157225daaf8SStefano Zampini 
1158225daaf8SStefano Zampini    Collective
1159225daaf8SStefano Zampini 
1160225daaf8SStefano Zampini    Input Parameters:
1161225daaf8SStefano Zampini +  vparcsr  - the pointer to the hypre_ParCSRMatrix
1162bb4689ddSStefano Zampini .  mtype    - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported.
1163225daaf8SStefano Zampini -  copymode - PETSc copying options
1164225daaf8SStefano Zampini 
1165225daaf8SStefano Zampini    Output Parameter:
1166225daaf8SStefano Zampini .  A  - the matrix
1167225daaf8SStefano Zampini 
1168225daaf8SStefano Zampini    Level: intermediate
1169225daaf8SStefano Zampini 
1170225daaf8SStefano Zampini .seealso: MatHYPRE, PetscCopyMode
1171225daaf8SStefano Zampini */
1172978814f1SStefano Zampini #undef __FUNCT__
1173225daaf8SStefano Zampini #define __FUNCT__ "MatCreateFromParCSR"
1174dd9c0a25Sstefano_zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(hypre_ParCSRMatrix *vparcsr, MatType mtype, PetscCopyMode copymode, Mat* A)
1175978814f1SStefano Zampini {
1176225daaf8SStefano Zampini   Mat                   T;
1177978814f1SStefano Zampini   Mat_HYPRE             *hA;
117863c07aadSStefano Zampini   hypre_ParCSRMatrix    *parcsr;
1179978814f1SStefano Zampini   MPI_Comm              comm;
1180978814f1SStefano Zampini   PetscInt              rstart,rend,cstart,cend,M,N;
1181bb4689ddSStefano Zampini   PetscBool             isseqaij,ismpiaij,isaij,ishyp,isis;
1182978814f1SStefano Zampini   PetscErrorCode        ierr;
1183978814f1SStefano Zampini 
1184978814f1SStefano Zampini   PetscFunctionBegin;
1185978814f1SStefano Zampini   parcsr = (hypre_ParCSRMatrix *)vparcsr;
1186978814f1SStefano Zampini   comm   = hypre_ParCSRMatrixComm(parcsr);
1187225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
1188225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
1189225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr);
1190225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr);
11918cfe8d00SStefano Zampini   ierr   = PetscStrcmp(mtype,MATIS,&isis);CHKERRQ(ierr);
1192225daaf8SStefano Zampini   isaij  = (PetscBool)(isseqaij || ismpiaij || isaij);
1193bb4689ddSStefano 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);
1194225daaf8SStefano Zampini   if (ishyp && copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES");
1195978814f1SStefano Zampini 
1196978814f1SStefano Zampini   /* access ParCSRMatrix */
1197978814f1SStefano Zampini   rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr);
1198978814f1SStefano Zampini   rend   = hypre_ParCSRMatrixLastRowIndex(parcsr);
1199978814f1SStefano Zampini   cstart = hypre_ParCSRMatrixFirstColDiag(parcsr);
1200978814f1SStefano Zampini   cend   = hypre_ParCSRMatrixLastColDiag(parcsr);
1201978814f1SStefano Zampini   M      = hypre_ParCSRMatrixGlobalNumRows(parcsr);
1202978814f1SStefano Zampini   N      = hypre_ParCSRMatrixGlobalNumCols(parcsr);
1203978814f1SStefano Zampini 
1204fa92c42cSstefano_zampini   /* fix for empty local rows/columns */
1205fa92c42cSstefano_zampini   if (rend < rstart) rend = rstart;
1206fa92c42cSstefano_zampini   if (cend < cstart) cend = cstart;
1207fa92c42cSstefano_zampini 
1208978814f1SStefano Zampini   /* create PETSc matrix with MatHYPRE */
1209225daaf8SStefano Zampini   ierr = MatCreate(comm,&T);CHKERRQ(ierr);
1210225daaf8SStefano Zampini   ierr = MatSetSizes(T,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr);
1211225daaf8SStefano Zampini   ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr);
1212225daaf8SStefano Zampini   hA   = (Mat_HYPRE*)(T->data);
1213978814f1SStefano Zampini 
1214978814f1SStefano Zampini   /* create HYPRE_IJMatrix */
1215978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
1216978814f1SStefano Zampini 
1217978814f1SStefano Zampini   /* set ParCSR object */
1218978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1219978814f1SStefano Zampini   hypre_IJMatrixObject(hA->ij) = parcsr;
12204ec6421dSstefano_zampini   T->preallocated = PETSC_TRUE;
1221978814f1SStefano Zampini 
1222978814f1SStefano Zampini   /* set assembled flag */
1223978814f1SStefano Zampini   hypre_IJMatrixAssembleFlag(hA->ij) = 1;
1224978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1225225daaf8SStefano Zampini   if (ishyp) {
12266d2a658fSstefano_zampini     PetscMPIInt myid = 0;
12276d2a658fSstefano_zampini 
12286d2a658fSstefano_zampini     /* make sure we always have row_starts and col_starts available */
12296d2a658fSstefano_zampini     if (HYPRE_AssumedPartitionCheck()) {
12306d2a658fSstefano_zampini       ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
12316d2a658fSstefano_zampini     }
12326d2a658fSstefano_zampini     if (!hypre_ParCSRMatrixOwnsColStarts(parcsr)) {
12336d2a658fSstefano_zampini       PetscLayout map;
12346d2a658fSstefano_zampini 
12356d2a658fSstefano_zampini       ierr = MatGetLayouts(T,NULL,&map);CHKERRQ(ierr);
12366d2a658fSstefano_zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
12376d2a658fSstefano_zampini       hypre_ParCSRMatrixColStarts(parcsr) = map->range + myid;
12386d2a658fSstefano_zampini     }
12396d2a658fSstefano_zampini     if (!hypre_ParCSRMatrixOwnsRowStarts(parcsr)) {
12406d2a658fSstefano_zampini       PetscLayout map;
12416d2a658fSstefano_zampini 
12426d2a658fSstefano_zampini       ierr = MatGetLayouts(T,&map,NULL);CHKERRQ(ierr);
12436d2a658fSstefano_zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
12446d2a658fSstefano_zampini       hypre_ParCSRMatrixRowStarts(parcsr) = map->range + myid;
12456d2a658fSstefano_zampini     }
1246978814f1SStefano Zampini     /* prevent from freeing the pointer */
1247978814f1SStefano Zampini     if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE;
1248225daaf8SStefano Zampini     *A   = T;
1249978814f1SStefano Zampini     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1250978814f1SStefano Zampini     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1251bb4689ddSStefano Zampini   } else if (isaij) {
1252bb4689ddSStefano Zampini     if (copymode != PETSC_OWN_POINTER) {
1253225daaf8SStefano Zampini       /* prevent from freeing the pointer */
1254225daaf8SStefano Zampini       hA->inner_free = PETSC_FALSE;
1255225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
1256225daaf8SStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
1257225daaf8SStefano Zampini     } else { /* AIJ return type with PETSC_OWN_POINTER */
1258225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1259225daaf8SStefano Zampini       *A   = T;
1260225daaf8SStefano Zampini     }
1261bb4689ddSStefano Zampini   } else if (isis) {
12628cfe8d00SStefano Zampini     ierr = MatConvert_HYPRE_IS(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
12638cfe8d00SStefano Zampini     if (copymode != PETSC_OWN_POINTER) hA->inner_free = PETSC_FALSE;
1264bb4689ddSStefano Zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
1265bb4689ddSStefano Zampini   }
1266978814f1SStefano Zampini   PetscFunctionReturn(0);
1267978814f1SStefano Zampini }
1268978814f1SStefano Zampini 
126963c07aadSStefano Zampini #undef __FUNCT__
1270dd9c0a25Sstefano_zampini #define __FUNCT__ "MatHYPREGetParCSR_HYPRE"
1271dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR_HYPRE(Mat A, hypre_ParCSRMatrix **parcsr)
1272dd9c0a25Sstefano_zampini {
1273dd9c0a25Sstefano_zampini   Mat_HYPRE*            hA = (Mat_HYPRE*)A->data;
1274dd9c0a25Sstefano_zampini   HYPRE_Int             type;
1275dd9c0a25Sstefano_zampini   PetscErrorCode        ierr;
1276dd9c0a25Sstefano_zampini 
1277dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1278dd9c0a25Sstefano_zampini   if (!hA->ij) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"HYPRE_IJMatrix not present");
1279dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
1280dd9c0a25Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"HYPRE_IJMatrix is not of type HYPRE_PARCSR");
1281dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)parcsr));
1282dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1283dd9c0a25Sstefano_zampini }
1284dd9c0a25Sstefano_zampini 
1285dd9c0a25Sstefano_zampini /*
1286dd9c0a25Sstefano_zampini    MatHYPREGetParCSR - Gets the pointer to the ParCSR matrix
1287dd9c0a25Sstefano_zampini 
1288dd9c0a25Sstefano_zampini    Not collective
1289dd9c0a25Sstefano_zampini 
1290dd9c0a25Sstefano_zampini    Input Parameters:
1291dd9c0a25Sstefano_zampini +  A  - the MATHYPRE object
1292dd9c0a25Sstefano_zampini 
1293dd9c0a25Sstefano_zampini    Output Parameter:
1294dd9c0a25Sstefano_zampini .  parcsr  - the pointer to the hypre_ParCSRMatrix
1295dd9c0a25Sstefano_zampini 
1296dd9c0a25Sstefano_zampini    Level: intermediate
1297dd9c0a25Sstefano_zampini 
1298dd9c0a25Sstefano_zampini .seealso: MatHYPRE, PetscCopyMode
1299dd9c0a25Sstefano_zampini */
1300dd9c0a25Sstefano_zampini #undef __FUNCT__
1301dd9c0a25Sstefano_zampini #define __FUNCT__ "MatHYPREGetParCSR"
1302dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR(Mat A, hypre_ParCSRMatrix **parcsr)
1303dd9c0a25Sstefano_zampini {
1304dd9c0a25Sstefano_zampini   PetscErrorCode ierr;
1305dd9c0a25Sstefano_zampini 
1306dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1307dd9c0a25Sstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1308dd9c0a25Sstefano_zampini   PetscValidType(A,1);
1309dd9c0a25Sstefano_zampini   ierr = PetscUseMethod(A,"MatHYPREGetParCSR_C",(Mat,hypre_ParCSRMatrix**),(A,parcsr));CHKERRQ(ierr);
1310dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1311dd9c0a25Sstefano_zampini }
1312dd9c0a25Sstefano_zampini 
1313dd9c0a25Sstefano_zampini #undef __FUNCT__
131463c07aadSStefano Zampini #define __FUNCT__ "MatCreate_HYPRE"
131563c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
131663c07aadSStefano Zampini {
131763c07aadSStefano Zampini   Mat_HYPRE      *hB;
131863c07aadSStefano Zampini   PetscErrorCode ierr;
131963c07aadSStefano Zampini 
132063c07aadSStefano Zampini   PetscFunctionBegin;
132163c07aadSStefano Zampini   ierr           = PetscNewLog(B,&hB);CHKERRQ(ierr);
1322978814f1SStefano Zampini   hB->inner_free = PETSC_TRUE;
1323978814f1SStefano Zampini 
132463c07aadSStefano Zampini   B->data       = (void*)hB;
132563c07aadSStefano Zampini   B->rmap->bs   = 1;
132663c07aadSStefano Zampini   B->assembled  = PETSC_FALSE;
132763c07aadSStefano Zampini 
132863c07aadSStefano Zampini   B->ops->mult          = MatMult_HYPRE;
132963c07aadSStefano Zampini   B->ops->multtranspose = MatMultTranspose_HYPRE;
133063c07aadSStefano Zampini   B->ops->setup         = MatSetUp_HYPRE;
133163c07aadSStefano Zampini   B->ops->destroy       = MatDestroy_HYPRE;
133263c07aadSStefano Zampini   B->ops->assemblyend   = MatAssemblyEnd_HYPRE;
1333cd8bc7baSStefano Zampini   B->ops->ptap          = MatPtAP_HYPRE_HYPRE;
1334d975228cSstefano_zampini   B->ops->setvalues     = MatSetValues_HYPRE;
133563c07aadSStefano Zampini 
133663c07aadSStefano Zampini   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
133763c07aadSStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
133863c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
13392df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr);
1340c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1341c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1342d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPRESetPreallocation_C",MatHYPRESetPreallocation_HYPRE);CHKERRQ(ierr);
1343dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPREGetParCSR_C",MatHYPREGetParCSR_HYPRE);CHKERRQ(ierr);
134463c07aadSStefano Zampini   PetscFunctionReturn(0);
134563c07aadSStefano Zampini }
134663c07aadSStefano Zampini 
1347225daaf8SStefano Zampini #undef __FUNCT__
1348225daaf8SStefano Zampini #define __FUNCT__ "hypre_array_destroy"
1349225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void *ptr)
1350225daaf8SStefano Zampini {
1351225daaf8SStefano Zampini    PetscFunctionBegin;
1352225daaf8SStefano Zampini    hypre_TFree(ptr);
1353225daaf8SStefano Zampini    PetscFunctionReturn(0);
1354225daaf8SStefano Zampini }
1355225daaf8SStefano Zampini 
135663c07aadSStefano Zampini #if 0
135763c07aadSStefano Zampini /*
135863c07aadSStefano Zampini     Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format
135963c07aadSStefano Zampini 
136063c07aadSStefano Zampini     This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first
136163c07aadSStefano Zampini     which will corrupt the PETSc data structure if we did this. Need a work around to this problem.
136263c07aadSStefano Zampini */
136363c07aadSStefano Zampini #include <_hypre_IJ_mv.h>
136463c07aadSStefano Zampini #include <HYPRE_IJ_mv.h>
136563c07aadSStefano Zampini 
136663c07aadSStefano Zampini #undef __FUNCT__
136763c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixLink"
136863c07aadSStefano Zampini PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij)
136963c07aadSStefano Zampini {
137063c07aadSStefano Zampini   PetscErrorCode        ierr;
137163c07aadSStefano Zampini   PetscInt              rstart,rend,cstart,cend;
137263c07aadSStefano Zampini   PetscBool             flg;
137363c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
137463c07aadSStefano Zampini 
137563c07aadSStefano Zampini   PetscFunctionBegin;
137663c07aadSStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
137763c07aadSStefano Zampini   PetscValidType(A,1);
137863c07aadSStefano Zampini   PetscValidPointer(ij,2);
137963c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
138063c07aadSStefano Zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices");
138163c07aadSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
138263c07aadSStefano Zampini 
138363c07aadSStefano Zampini   ierr   = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
138463c07aadSStefano Zampini   rstart = A->rmap->rstart;
138563c07aadSStefano Zampini   rend   = A->rmap->rend;
138663c07aadSStefano Zampini   cstart = A->cmap->rstart;
138763c07aadSStefano Zampini   cend   = A->cmap->rend;
138863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij));
138963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
139063c07aadSStefano Zampini 
139163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij));
139263c07aadSStefano Zampini   PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij));
139363c07aadSStefano Zampini 
139463c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
139563c07aadSStefano Zampini 
139663c07aadSStefano Zampini   /* this is the Hack part where we monkey directly with the hypre datastructures */
139763c07aadSStefano Zampini 
139863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij));
139963c07aadSStefano Zampini   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
140063c07aadSStefano Zampini   PetscFunctionReturn(0);
140163c07aadSStefano Zampini }
140263c07aadSStefano Zampini #endif
1403