xref: /petsc/src/mat/impls/hypre/mhypre.c (revision a055b5aa40477f4e4106d726f8244c01e023706a)
163c07aadSStefano Zampini 
263c07aadSStefano Zampini /*
363c07aadSStefano Zampini     Creates hypre ijmatrix from PETSc matrix
463c07aadSStefano Zampini */
5225daaf8SStefano Zampini 
6dd9c0a25Sstefano_zampini #include <petscmathypre.h>
763c07aadSStefano Zampini #include <petsc/private/matimpl.h>
863c07aadSStefano Zampini #include <../src/mat/impls/hypre/mhypre.h>
963c07aadSStefano Zampini #include <../src/mat/impls/aij/mpi/mpiaij.h>
1058968eb6SStefano Zampini #include <../src/vec/vec/impls/hypre/vhyp.h>
1158968eb6SStefano Zampini #include <HYPRE.h>
12c1a070e6SStefano Zampini #include <HYPRE_utilities.h>
13cd8bc7baSStefano Zampini #include <_hypre_parcsr_ls.h>
1463c07aadSStefano Zampini 
1563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*);
1663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix);
1763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
1863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
1963c07aadSStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool);
20225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void*);
2163c07aadSStefano Zampini 
2263c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij)
2363c07aadSStefano Zampini {
2463c07aadSStefano Zampini   PetscErrorCode ierr;
2563c07aadSStefano Zampini   PetscInt       i,n_d,n_o;
2663c07aadSStefano Zampini   const PetscInt *ia_d,*ia_o;
2763c07aadSStefano Zampini   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
2863c07aadSStefano Zampini   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
2963c07aadSStefano Zampini 
3063c07aadSStefano Zampini   PetscFunctionBegin;
3163c07aadSStefano Zampini   if (A_d) { /* determine number of nonzero entries in local diagonal part */
3263c07aadSStefano Zampini     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
3363c07aadSStefano Zampini     if (done_d) {
3463c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr);
3563c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
3663c07aadSStefano Zampini         nnz_d[i] = ia_d[i+1] - ia_d[i];
3763c07aadSStefano Zampini       }
3863c07aadSStefano Zampini     }
3963c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
4063c07aadSStefano Zampini   }
4163c07aadSStefano Zampini   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
4263c07aadSStefano Zampini     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
4363c07aadSStefano Zampini     if (done_o) {
4463c07aadSStefano Zampini       ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr);
4563c07aadSStefano Zampini       for (i=0; i<n_o; i++) {
4663c07aadSStefano Zampini         nnz_o[i] = ia_o[i+1] - ia_o[i];
4763c07aadSStefano Zampini       }
4863c07aadSStefano Zampini     }
4963c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
5063c07aadSStefano Zampini   }
5163c07aadSStefano Zampini   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
5263c07aadSStefano Zampini     if (!done_o) { /* only diagonal part */
5363c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
5463c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
5563c07aadSStefano Zampini         nnz_o[i] = 0;
5663c07aadSStefano Zampini       }
5763c07aadSStefano Zampini     }
5863c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o));
5963c07aadSStefano Zampini     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
6063c07aadSStefano Zampini     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
6163c07aadSStefano Zampini   }
6263c07aadSStefano Zampini   PetscFunctionReturn(0);
6363c07aadSStefano Zampini }
6463c07aadSStefano Zampini 
6563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA)
6663c07aadSStefano Zampini {
6763c07aadSStefano Zampini   PetscErrorCode ierr;
6863c07aadSStefano Zampini   PetscInt       rstart,rend,cstart,cend;
6963c07aadSStefano Zampini 
7063c07aadSStefano Zampini   PetscFunctionBegin;
71d975228cSstefano_zampini   ierr   = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
72d975228cSstefano_zampini   ierr   = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
7363c07aadSStefano Zampini   rstart = A->rmap->rstart;
7463c07aadSStefano Zampini   rend   = A->rmap->rend;
7563c07aadSStefano Zampini   cstart = A->cmap->rstart;
7663c07aadSStefano Zampini   cend   = A->cmap->rend;
7763c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
7863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
7963c07aadSStefano Zampini   {
8063c07aadSStefano Zampini     PetscBool      same;
8163c07aadSStefano Zampini     Mat            A_d,A_o;
8263c07aadSStefano Zampini     const PetscInt *colmap;
8363c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
8463c07aadSStefano Zampini     if (same) {
8563c07aadSStefano Zampini       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
8663c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
8763c07aadSStefano Zampini       PetscFunctionReturn(0);
8863c07aadSStefano Zampini     }
8963c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
9063c07aadSStefano Zampini     if (same) {
9163c07aadSStefano Zampini       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
9263c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
9363c07aadSStefano Zampini       PetscFunctionReturn(0);
9463c07aadSStefano Zampini     }
9563c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
9663c07aadSStefano Zampini     if (same) {
9763c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
9863c07aadSStefano Zampini       PetscFunctionReturn(0);
9963c07aadSStefano Zampini     }
10063c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
10163c07aadSStefano Zampini     if (same) {
10263c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
10363c07aadSStefano Zampini       PetscFunctionReturn(0);
10463c07aadSStefano Zampini     }
10563c07aadSStefano Zampini   }
10663c07aadSStefano Zampini   PetscFunctionReturn(0);
10763c07aadSStefano Zampini }
10863c07aadSStefano Zampini 
10963c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij)
11063c07aadSStefano Zampini {
11163c07aadSStefano Zampini   PetscErrorCode    ierr;
11263c07aadSStefano Zampini   PetscInt          i,rstart,rend,ncols,nr,nc;
11363c07aadSStefano Zampini   const PetscScalar *values;
11463c07aadSStefano Zampini   const PetscInt    *cols;
11563c07aadSStefano Zampini   PetscBool         flg;
11663c07aadSStefano Zampini 
11763c07aadSStefano Zampini   PetscFunctionBegin;
11863c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
11963c07aadSStefano Zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
12063c07aadSStefano Zampini   if (flg && nr == nc) {
12163c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
12263c07aadSStefano Zampini     PetscFunctionReturn(0);
12363c07aadSStefano Zampini   }
12463c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
12563c07aadSStefano Zampini   if (flg) {
12663c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
12763c07aadSStefano Zampini     PetscFunctionReturn(0);
12863c07aadSStefano Zampini   }
12963c07aadSStefano Zampini 
13063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
13163c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
13263c07aadSStefano Zampini   for (i=rstart; i<rend; i++) {
13363c07aadSStefano Zampini     ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
134e3977e59Sstefano_zampini     if (ncols) {
13563c07aadSStefano Zampini       PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
136e3977e59Sstefano_zampini     }
13763c07aadSStefano Zampini     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
13863c07aadSStefano Zampini   }
13963c07aadSStefano Zampini   PetscFunctionReturn(0);
14063c07aadSStefano Zampini }
14163c07aadSStefano Zampini 
14263c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij)
14363c07aadSStefano Zampini {
14463c07aadSStefano Zampini   PetscErrorCode        ierr;
14563c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag = (Mat_SeqAIJ*)A->data;
14658968eb6SStefano Zampini   HYPRE_Int             type;
14763c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
14863c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
14963c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag;
15063c07aadSStefano Zampini 
15163c07aadSStefano Zampini   PetscFunctionBegin;
15263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
153ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
154ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
155ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
15663c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
15763c07aadSStefano Zampini   /*
15863c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
15963c07aadSStefano Zampini   */
16063c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));
16163c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));
16263c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
163ea9daf28SStefano Zampini 
164ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
16563c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
16663c07aadSStefano Zampini   PetscFunctionReturn(0);
16763c07aadSStefano Zampini }
16863c07aadSStefano Zampini 
16963c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
17063c07aadSStefano Zampini {
17163c07aadSStefano Zampini   PetscErrorCode        ierr;
17263c07aadSStefano Zampini   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
17363c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag,*poffd;
17463c07aadSStefano Zampini   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
17558968eb6SStefano Zampini   HYPRE_Int             type;
17663c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
17763c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
17863c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag,*hoffd;
17963c07aadSStefano Zampini 
18063c07aadSStefano Zampini   PetscFunctionBegin;
18163c07aadSStefano Zampini   pdiag = (Mat_SeqAIJ*) pA->A->data;
18263c07aadSStefano Zampini   poffd = (Mat_SeqAIJ*) pA->B->data;
18363c07aadSStefano Zampini   /* cstart is only valid for square MPIAIJ layed out in the usual way */
18463c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
18563c07aadSStefano Zampini 
18663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
187ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
188ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
189ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
19063c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
19163c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
19263c07aadSStefano Zampini 
19363c07aadSStefano Zampini   /*
19463c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
19563c07aadSStefano Zampini   */
19663c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
19763c07aadSStefano Zampini   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
19863c07aadSStefano Zampini   jj  = (PetscInt*)hdiag->j;
19963c07aadSStefano Zampini   pjj = (PetscInt*)pdiag->j;
20063c07aadSStefano Zampini   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
20163c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
20263c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
20363c07aadSStefano Zampini   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
20463c07aadSStefano Zampini      If we hacked a hypre a bit more we might be able to avoid this step */
20563c07aadSStefano Zampini   jj  = (PetscInt*) hoffd->j;
20663c07aadSStefano Zampini   pjj = (PetscInt*) poffd->j;
20763c07aadSStefano Zampini   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
20863c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
20963c07aadSStefano Zampini 
210ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
21163c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
21263c07aadSStefano Zampini   PetscFunctionReturn(0);
21363c07aadSStefano Zampini }
21463c07aadSStefano Zampini 
2152df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B)
2162df22349SStefano Zampini {
2172df22349SStefano Zampini   Mat_HYPRE*             mhA = (Mat_HYPRE*)(A->data);
2182df22349SStefano Zampini   Mat                    lA;
2192df22349SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
2202df22349SStefano Zampini   IS                     is;
2212df22349SStefano Zampini   hypre_ParCSRMatrix     *hA;
2222df22349SStefano Zampini   hypre_CSRMatrix        *hdiag,*hoffd;
2232df22349SStefano Zampini   MPI_Comm               comm;
2242df22349SStefano Zampini   PetscScalar            *hdd,*hod,*aa,*data;
2257d968826Sstefano_zampini   HYPRE_Int              *col_map_offd,*hdi,*hdj,*hoi,*hoj;
2262df22349SStefano Zampini   PetscInt               *ii,*jj,*iptr,*jptr;
2272df22349SStefano Zampini   PetscInt               cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N;
22858968eb6SStefano Zampini   HYPRE_Int              type;
2292df22349SStefano Zampini   PetscErrorCode         ierr;
2302df22349SStefano Zampini 
2312df22349SStefano Zampini   PetscFunctionBegin;
232a1787963SStefano Zampini   comm = PetscObjectComm((PetscObject)A);
2332df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type));
2342df22349SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
2352df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA));
2362df22349SStefano Zampini   M     = hypre_ParCSRMatrixGlobalNumRows(hA);
2372df22349SStefano Zampini   N     = hypre_ParCSRMatrixGlobalNumCols(hA);
2382df22349SStefano Zampini   str   = hypre_ParCSRMatrixFirstRowIndex(hA);
2392df22349SStefano Zampini   stc   = hypre_ParCSRMatrixFirstColDiag(hA);
2402df22349SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(hA);
2412df22349SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(hA);
2422df22349SStefano Zampini   dr    = hypre_CSRMatrixNumRows(hdiag);
2432df22349SStefano Zampini   dc    = hypre_CSRMatrixNumCols(hdiag);
2442df22349SStefano Zampini   nnz   = hypre_CSRMatrixNumNonzeros(hdiag);
2452df22349SStefano Zampini   hdi   = hypre_CSRMatrixI(hdiag);
2462df22349SStefano Zampini   hdj   = hypre_CSRMatrixJ(hdiag);
2472df22349SStefano Zampini   hdd   = hypre_CSRMatrixData(hdiag);
2482df22349SStefano Zampini   oc    = hypre_CSRMatrixNumCols(hoffd);
2492df22349SStefano Zampini   nnz  += hypre_CSRMatrixNumNonzeros(hoffd);
2502df22349SStefano Zampini   hoi   = hypre_CSRMatrixI(hoffd);
2512df22349SStefano Zampini   hoj   = hypre_CSRMatrixJ(hoffd);
2522df22349SStefano Zampini   hod   = hypre_CSRMatrixData(hoffd);
2532df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
2542df22349SStefano Zampini     PetscInt *aux;
2552df22349SStefano Zampini 
2562df22349SStefano Zampini     /* generate l2g maps for rows and cols */
2572df22349SStefano Zampini     ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr);
2582df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
2592df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2602df22349SStefano Zampini     col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
2612df22349SStefano Zampini     ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr);
2622df22349SStefano Zampini     for (i=0; i<dc; i++) aux[i] = i+stc;
2632df22349SStefano Zampini     for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i];
2642df22349SStefano Zampini     ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
2652df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
2662df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2672df22349SStefano Zampini     /* create MATIS object */
2682df22349SStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
2692df22349SStefano Zampini     ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr);
2702df22349SStefano Zampini     ierr = MatSetType(*B,MATIS);CHKERRQ(ierr);
2712df22349SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr);
2722df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
2732df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
2742df22349SStefano Zampini 
2752df22349SStefano Zampini     /* allocate CSR for local matrix */
2762df22349SStefano Zampini     ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr);
2772df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr);
2782df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr);
2792df22349SStefano Zampini   } else {
2802df22349SStefano Zampini     PetscInt  nr;
2812df22349SStefano Zampini     PetscBool done;
2822df22349SStefano Zampini     ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr);
2832df22349SStefano Zampini     ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr);
2842df22349SStefano 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);
2852df22349SStefano 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);
2862df22349SStefano Zampini     ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr);
2872df22349SStefano Zampini   }
2882df22349SStefano Zampini   /* merge local matrices */
2892df22349SStefano Zampini   ii   = iptr;
2902df22349SStefano Zampini   jj   = jptr;
2912df22349SStefano Zampini   aa   = data;
2922df22349SStefano Zampini   *ii  = *(hdi++) + *(hoi++);
2932df22349SStefano Zampini   for (jd=0,jo=0,cum=0; *ii<nnz; cum++) {
2942df22349SStefano Zampini     PetscScalar *aold = aa;
2952df22349SStefano Zampini     PetscInt    *jold = jj,nc = jd+jo;
2962df22349SStefano Zampini     for (; jd<*hdi; jd++) { *jj++ = *hdj++;      *aa++ = *hdd++; }
2972df22349SStefano Zampini     for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
2982df22349SStefano Zampini     *(++ii) = *(hdi++) + *(hoi++);
2992df22349SStefano Zampini     ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr);
3002df22349SStefano Zampini   }
3012df22349SStefano Zampini   for (; cum<dr; cum++) *(++ii) = nnz;
3022df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
303a033916dSStefano Zampini     Mat_SeqAIJ* a;
304a033916dSStefano Zampini 
3052df22349SStefano Zampini     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr);
3062df22349SStefano Zampini     ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr);
307a033916dSStefano Zampini     /* hack SeqAIJ */
308a033916dSStefano Zampini     a          = (Mat_SeqAIJ*)(lA->data);
309a033916dSStefano Zampini     a->free_a  = PETSC_TRUE;
310a033916dSStefano Zampini     a->free_ij = PETSC_TRUE;
3112df22349SStefano Zampini     ierr = MatDestroy(&lA);CHKERRQ(ierr);
3122df22349SStefano Zampini   }
3132df22349SStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3142df22349SStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3152df22349SStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
3162df22349SStefano Zampini     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
3172df22349SStefano Zampini   }
3182df22349SStefano Zampini   PetscFunctionReturn(0);
3192df22349SStefano Zampini }
3202df22349SStefano Zampini 
32163c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
32263c07aadSStefano Zampini {
32363c07aadSStefano Zampini   Mat_HYPRE      *hB;
32463c07aadSStefano Zampini   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
32563c07aadSStefano Zampini   PetscErrorCode ierr;
32663c07aadSStefano Zampini 
32763c07aadSStefano Zampini   PetscFunctionBegin;
32863c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
32963c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
33063c07aadSStefano Zampini     /* always destroy the old matrix and create a new memory;
33163c07aadSStefano Zampini        hope this does not churn the memory too much. The problem
33263c07aadSStefano Zampini        is I do not know if it is possible to put the matrix back to
33363c07aadSStefano Zampini        its initial state so that we can directly copy the values
33463c07aadSStefano Zampini        the second time through. */
33563c07aadSStefano Zampini     hB = (Mat_HYPRE*)((*B)->data);
33663c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
33763c07aadSStefano Zampini   } else {
33863c07aadSStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
33963c07aadSStefano Zampini     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
34063c07aadSStefano Zampini     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
34163c07aadSStefano Zampini     hB   = (Mat_HYPRE*)((*B)->data);
34263c07aadSStefano Zampini   }
34363c07aadSStefano Zampini   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
34463c07aadSStefano Zampini   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
3454ec6421dSstefano_zampini   (*B)->preallocated = PETSC_TRUE;
34663c07aadSStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
34763c07aadSStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
34863c07aadSStefano Zampini   PetscFunctionReturn(0);
34963c07aadSStefano Zampini }
35063c07aadSStefano Zampini 
351ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
35263c07aadSStefano Zampini {
35363c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
35463c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
35563c07aadSStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
35663c07aadSStefano Zampini   MPI_Comm           comm;
35763c07aadSStefano Zampini   PetscScalar        *da,*oa,*aptr;
35863c07aadSStefano Zampini   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
35963c07aadSStefano Zampini   PetscInt           i,dnnz,onnz,m,n;
36058968eb6SStefano Zampini   HYPRE_Int          type;
36163c07aadSStefano Zampini   PetscMPIInt        size;
36263c07aadSStefano Zampini   PetscErrorCode     ierr;
36363c07aadSStefano Zampini 
36463c07aadSStefano Zampini   PetscFunctionBegin;
36563c07aadSStefano Zampini   comm = PetscObjectComm((PetscObject)A);
36663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
36763c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
36863c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
36963c07aadSStefano Zampini     PetscBool ismpiaij,isseqaij;
37063c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
37163c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
37263c07aadSStefano Zampini     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
37363c07aadSStefano Zampini   }
37463c07aadSStefano Zampini   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
37563c07aadSStefano Zampini 
37663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
37763c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(parcsr);
37863c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(parcsr);
37963c07aadSStefano Zampini   m     = hypre_CSRMatrixNumRows(hdiag);
38063c07aadSStefano Zampini   n     = hypre_CSRMatrixNumCols(hdiag);
38163c07aadSStefano Zampini   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
38263c07aadSStefano Zampini   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
383225daaf8SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
38463c07aadSStefano Zampini     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
38563c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
38663c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
387225daaf8SStefano Zampini   } else if (reuse == MAT_REUSE_MATRIX) {
38863c07aadSStefano Zampini     PetscInt  nr;
38963c07aadSStefano Zampini     PetscBool done;
39063c07aadSStefano Zampini     if (size > 1) {
39163c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
39263c07aadSStefano Zampini 
39363c07aadSStefano Zampini       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
39463c07aadSStefano 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);
39563c07aadSStefano 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);
39663c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
39763c07aadSStefano Zampini     } else {
39863c07aadSStefano Zampini       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
39963c07aadSStefano Zampini       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
40063c07aadSStefano 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);
40163c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
40263c07aadSStefano Zampini     }
403225daaf8SStefano Zampini   } else { /* MAT_INPLACE_MATRIX */
4047d968826Sstefano_zampini     dii = (PetscInt*)hypre_CSRMatrixI(hdiag);
4057d968826Sstefano_zampini     djj = (PetscInt*)hypre_CSRMatrixJ(hdiag);
406225daaf8SStefano Zampini     da  = hypre_CSRMatrixData(hdiag);
40763c07aadSStefano Zampini   }
40863c07aadSStefano Zampini   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
40963c07aadSStefano Zampini   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
41063c07aadSStefano Zampini   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
41163c07aadSStefano Zampini   iptr = djj;
41263c07aadSStefano Zampini   aptr = da;
41363c07aadSStefano Zampini   for (i=0; i<m; i++) {
41463c07aadSStefano Zampini     PetscInt nc = dii[i+1]-dii[i];
41563c07aadSStefano Zampini     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
41663c07aadSStefano Zampini     iptr += nc;
41763c07aadSStefano Zampini     aptr += nc;
41863c07aadSStefano Zampini   }
41963c07aadSStefano Zampini   if (size > 1) {
4207d968826Sstefano_zampini     HYPRE_Int *offdj,*coffd;
42163c07aadSStefano Zampini 
422225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
42363c07aadSStefano Zampini       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
42463c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
42563c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
426225daaf8SStefano Zampini     } else if (reuse == MAT_REUSE_MATRIX) {
42763c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
42863c07aadSStefano Zampini       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
42963c07aadSStefano Zampini       PetscBool  done;
43063c07aadSStefano Zampini 
43163c07aadSStefano Zampini       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
43263c07aadSStefano 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);
43363c07aadSStefano 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);
43463c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
435225daaf8SStefano Zampini     } else { /* MAT_INPLACE_MATRIX */
4367d968826Sstefano_zampini       oii = (PetscInt*)hypre_CSRMatrixI(hoffd);
4377d968826Sstefano_zampini       ojj = (PetscInt*)hypre_CSRMatrixJ(hoffd);
438225daaf8SStefano Zampini       oa  = hypre_CSRMatrixData(hoffd);
43963c07aadSStefano Zampini     }
44063c07aadSStefano Zampini     ierr  = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
44163c07aadSStefano Zampini     offdj = hypre_CSRMatrixJ(hoffd);
44263c07aadSStefano Zampini     coffd = hypre_ParCSRMatrixColMapOffd(parcsr);
44363c07aadSStefano Zampini     for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]];
44463c07aadSStefano Zampini     ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr);
44563c07aadSStefano Zampini     iptr = ojj;
44663c07aadSStefano Zampini     aptr = oa;
44763c07aadSStefano Zampini     for (i=0; i<m; i++) {
44863c07aadSStefano Zampini        PetscInt nc = oii[i+1]-oii[i];
44963c07aadSStefano Zampini        ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
45063c07aadSStefano Zampini        iptr += nc;
45163c07aadSStefano Zampini        aptr += nc;
45263c07aadSStefano Zampini     }
453225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
45463c07aadSStefano Zampini       Mat_MPIAIJ *b;
45563c07aadSStefano Zampini       Mat_SeqAIJ *d,*o;
456225daaf8SStefano Zampini 
45741073d71Sstefano_zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,B);CHKERRQ(ierr);
45863c07aadSStefano Zampini       /* hack MPIAIJ */
45963c07aadSStefano Zampini       b          = (Mat_MPIAIJ*)((*B)->data);
46063c07aadSStefano Zampini       d          = (Mat_SeqAIJ*)b->A->data;
46163c07aadSStefano Zampini       o          = (Mat_SeqAIJ*)b->B->data;
46263c07aadSStefano Zampini       d->free_a  = PETSC_TRUE;
46363c07aadSStefano Zampini       d->free_ij = PETSC_TRUE;
46463c07aadSStefano Zampini       o->free_a  = PETSC_TRUE;
46563c07aadSStefano Zampini       o->free_ij = PETSC_TRUE;
466225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
467225daaf8SStefano Zampini       Mat T;
46841073d71Sstefano_zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,&T);CHKERRQ(ierr);
469225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
470225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
471225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
472225daaf8SStefano Zampini       hypre_CSRMatrixI(hoffd)    = NULL;
473225daaf8SStefano Zampini       hypre_CSRMatrixJ(hoffd)    = NULL;
474225daaf8SStefano Zampini       hypre_CSRMatrixData(hoffd) = NULL;
475225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
47663c07aadSStefano Zampini     }
477225daaf8SStefano Zampini   } else {
478225daaf8SStefano Zampini     oii  = NULL;
479225daaf8SStefano Zampini     ojj  = NULL;
480225daaf8SStefano Zampini     oa   = NULL;
481225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
48263c07aadSStefano Zampini       Mat_SeqAIJ* b;
48363c07aadSStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
48463c07aadSStefano Zampini       /* hack SeqAIJ */
48563c07aadSStefano Zampini       b          = (Mat_SeqAIJ*)((*B)->data);
48663c07aadSStefano Zampini       b->free_a  = PETSC_TRUE;
48763c07aadSStefano Zampini       b->free_ij = PETSC_TRUE;
488225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
489225daaf8SStefano Zampini       Mat T;
490225daaf8SStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,&T);CHKERRQ(ierr);
491225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
492225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
493225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
494225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
49563c07aadSStefano Zampini     }
496225daaf8SStefano Zampini   }
497225daaf8SStefano Zampini 
498225daaf8SStefano Zampini   /* we have to use hypre_Tfree to free the arrays */
49963c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
500225daaf8SStefano Zampini     void *ptrs[6] = {dii,djj,da,oii,ojj,oa};
501225daaf8SStefano Zampini     const char *names[6] = {"_hypre_csr_dii",
502225daaf8SStefano Zampini                             "_hypre_csr_djj",
503225daaf8SStefano Zampini                             "_hypre_csr_da",
504225daaf8SStefano Zampini                             "_hypre_csr_oii",
505225daaf8SStefano Zampini                             "_hypre_csr_ojj",
506225daaf8SStefano Zampini                             "_hypre_csr_oa"};
507225daaf8SStefano Zampini     for (i=0; i<6; i++) {
508225daaf8SStefano Zampini       PetscContainer c;
509225daaf8SStefano Zampini 
510225daaf8SStefano Zampini       ierr = PetscContainerCreate(comm,&c);CHKERRQ(ierr);
511225daaf8SStefano Zampini       ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr);
512225daaf8SStefano Zampini       ierr = PetscContainerSetUserDestroy(c,hypre_array_destroy);CHKERRQ(ierr);
513225daaf8SStefano Zampini       ierr = PetscObjectCompose((PetscObject)(*B),names[i],(PetscObject)c);CHKERRQ(ierr);
514225daaf8SStefano Zampini       ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
515225daaf8SStefano Zampini     }
51663c07aadSStefano Zampini   }
51763c07aadSStefano Zampini   PetscFunctionReturn(0);
51863c07aadSStefano Zampini }
51963c07aadSStefano Zampini 
520613e5ff0Sstefano_zampini static PetscErrorCode MatAIJGetParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
521c1a070e6SStefano Zampini {
522613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *tA;
523c1a070e6SStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
524c1a070e6SStefano Zampini   Mat_SeqAIJ         *diag,*offd;
525c1a070e6SStefano Zampini   PetscInt           *garray,noffd,dnnz,onnz,*row_starts,*col_starts;
526c1a070e6SStefano Zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
527613e5ff0Sstefano_zampini   PetscBool          ismpiaij,isseqaij;
528c1a070e6SStefano Zampini   PetscErrorCode     ierr;
529c1a070e6SStefano Zampini 
530c1a070e6SStefano Zampini   PetscFunctionBegin;
531c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
532c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
533c1a070e6SStefano Zampini   if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type);
534c1a070e6SStefano Zampini   if (ismpiaij) {
535c1a070e6SStefano Zampini     Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data);
536c1a070e6SStefano Zampini 
537c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)a->A->data;
538c1a070e6SStefano Zampini     offd   = (Mat_SeqAIJ*)a->B->data;
539c1a070e6SStefano Zampini     garray = a->garray;
540c1a070e6SStefano Zampini     noffd  = a->B->cmap->N;
541c1a070e6SStefano Zampini     dnnz   = diag->nz;
542c1a070e6SStefano Zampini     onnz   = offd->nz;
543c1a070e6SStefano Zampini   } else {
544c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)A->data;
545c1a070e6SStefano Zampini     offd   = NULL;
546c1a070e6SStefano Zampini     garray = NULL;
547c1a070e6SStefano Zampini     noffd  = 0;
548c1a070e6SStefano Zampini     dnnz   = diag->nz;
549c1a070e6SStefano Zampini     onnz   = 0;
550c1a070e6SStefano Zampini   }
551225daaf8SStefano Zampini 
552c1a070e6SStefano Zampini   /* create a temporary ParCSR */
553c1a070e6SStefano Zampini   if (HYPRE_AssumedPartitionCheck()) {
554c1a070e6SStefano Zampini     PetscMPIInt myid;
555c1a070e6SStefano Zampini 
556c1a070e6SStefano Zampini     ierr       = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
557c1a070e6SStefano Zampini     row_starts = A->rmap->range + myid;
558c1a070e6SStefano Zampini     col_starts = A->cmap->range + myid;
559c1a070e6SStefano Zampini   } else {
560c1a070e6SStefano Zampini     row_starts = A->rmap->range;
561c1a070e6SStefano Zampini     col_starts = A->cmap->range;
562c1a070e6SStefano Zampini   }
5637d968826Sstefano_zampini   tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,(HYPRE_Int*)row_starts,(HYPRE_Int*)col_starts,noffd,dnnz,onnz);
564c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(tA,0);
565c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(tA,0);
566c1a070e6SStefano Zampini 
567225daaf8SStefano Zampini   /* set diagonal part */
568c1a070e6SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(tA);
5697d968826Sstefano_zampini   hypre_CSRMatrixI(hdiag)           = (HYPRE_Int*)diag->i;
5707d968826Sstefano_zampini   hypre_CSRMatrixJ(hdiag)           = (HYPRE_Int*)diag->j;
571c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = diag->a;
572c1a070e6SStefano Zampini   hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz;
573c1a070e6SStefano Zampini   hypre_CSRMatrixSetRownnz(hdiag);
574c1a070e6SStefano Zampini   hypre_CSRMatrixSetDataOwner(hdiag,0);
575c1a070e6SStefano Zampini 
576225daaf8SStefano Zampini   /* set offdiagonal part */
577c1a070e6SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(tA);
578c1a070e6SStefano Zampini   if (offd) {
5797d968826Sstefano_zampini     hypre_CSRMatrixI(hoffd)           = (HYPRE_Int*)offd->i;
5807d968826Sstefano_zampini     hypre_CSRMatrixJ(hoffd)           = (HYPRE_Int*)offd->j;
581c1a070e6SStefano Zampini     hypre_CSRMatrixData(hoffd)        = offd->a;
582c1a070e6SStefano Zampini     hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz;
583c1a070e6SStefano Zampini     hypre_CSRMatrixSetRownnz(hoffd);
584c1a070e6SStefano Zampini     hypre_CSRMatrixSetDataOwner(hoffd,0);
585c1a070e6SStefano Zampini     hypre_ParCSRMatrixSetNumNonzeros(tA);
5867d968826Sstefano_zampini     hypre_ParCSRMatrixColMapOffd(tA) = (HYPRE_Int*)garray;
587c1a070e6SStefano Zampini   }
588613e5ff0Sstefano_zampini   *hA = tA;
589613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
590613e5ff0Sstefano_zampini }
591c1a070e6SStefano Zampini 
592613e5ff0Sstefano_zampini static PetscErrorCode MatAIJRestoreParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
593613e5ff0Sstefano_zampini {
594613e5ff0Sstefano_zampini   hypre_CSRMatrix    *hdiag,*hoffd;
595c1a070e6SStefano Zampini 
596613e5ff0Sstefano_zampini   PetscFunctionBegin;
597613e5ff0Sstefano_zampini   hdiag = hypre_ParCSRMatrixDiag(*hA);
598613e5ff0Sstefano_zampini   hoffd = hypre_ParCSRMatrixOffd(*hA);
599c1a070e6SStefano Zampini   /* set pointers to NULL before destroying tA */
600c1a070e6SStefano Zampini   hypre_CSRMatrixI(hdiag)           = NULL;
601c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hdiag)           = NULL;
602c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = NULL;
603c1a070e6SStefano Zampini   hypre_CSRMatrixI(hoffd)           = NULL;
604c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hoffd)           = NULL;
605c1a070e6SStefano Zampini   hypre_CSRMatrixData(hoffd)        = NULL;
606613e5ff0Sstefano_zampini   hypre_ParCSRMatrixColMapOffd(*hA) = NULL;
607613e5ff0Sstefano_zampini   hypre_ParCSRMatrixDestroy(*hA);
608613e5ff0Sstefano_zampini   *hA = NULL;
609613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
610613e5ff0Sstefano_zampini }
611613e5ff0Sstefano_zampini 
612613e5ff0Sstefano_zampini /* calls RAP from BoomerAMG:
6133dad0653Sstefano_zampini    the resulting ParCSR will not own the column and row starts
6143dad0653Sstefano_zampini    It looks like we don't need to have the diagonal entries
6153dad0653Sstefano_zampini    ordered first in the rows of the diagonal part
6163dad0653Sstefano_zampini    for boomerAMGBuildCoarseOperator to work */
617*a055b5aaSBarry Smith static PetscErrorCode MatHYPRE_ParCSR_RAP(hypre_ParCSRMatrix *hR, hypre_ParCSRMatrix *hA,hypre_ParCSRMatrix *hP, hypre_ParCSRMatrix **hRAP)
618613e5ff0Sstefano_zampini {
619613e5ff0Sstefano_zampini   PetscErrorCode ierr;
620613e5ff0Sstefano_zampini   HYPRE_Int      P_owns_col_starts,R_owns_row_starts;
621613e5ff0Sstefano_zampini 
622613e5ff0Sstefano_zampini   PetscFunctionBegin;
623613e5ff0Sstefano_zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(hP);
624613e5ff0Sstefano_zampini   R_owns_row_starts = hypre_ParCSRMatrixOwnsRowStarts(hR);
625613e5ff0Sstefano_zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(hR,hA,hP,hRAP));
626613e5ff0Sstefano_zampini   PetscStackCallStandard(hypre_ParCSRMatrixSetNumNonzeros,(*hRAP));
627613e5ff0Sstefano_zampini   /* hypre_BoomerAMGBuildCoarseOperator steals the col_starts from P and the row_starts from R */
628613e5ff0Sstefano_zampini   hypre_ParCSRMatrixSetRowStartsOwner(*hRAP,0);
629613e5ff0Sstefano_zampini   hypre_ParCSRMatrixSetColStartsOwner(*hRAP,0);
630613e5ff0Sstefano_zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(hP,1);
631613e5ff0Sstefano_zampini   if (R_owns_row_starts) hypre_ParCSRMatrixSetRowStartsOwner(hR,1);
632613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
633613e5ff0Sstefano_zampini }
634613e5ff0Sstefano_zampini 
6356f231fbdSstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat P,Mat C)
636613e5ff0Sstefano_zampini {
6376f231fbdSstefano_zampini   Mat                B;
638613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *hA,*hP,*hPtAP;
639613e5ff0Sstefano_zampini   PetscErrorCode     ierr;
640613e5ff0Sstefano_zampini 
641613e5ff0Sstefano_zampini   PetscFunctionBegin;
642613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
643613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(P,&hP);CHKERRQ(ierr);
644613e5ff0Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(hP,hA,hP,&hPtAP);CHKERRQ(ierr);
6456f231fbdSstefano_zampini   ierr = MatCreateFromParCSR(hPtAP,MATAIJ,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
6466f231fbdSstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
647613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
648613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(P,&hP);CHKERRQ(ierr);
6496f231fbdSstefano_zampini   PetscFunctionReturn(0);
6506f231fbdSstefano_zampini }
6516f231fbdSstefano_zampini 
6526f231fbdSstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat P,PetscReal fill,Mat *C)
6536f231fbdSstefano_zampini {
6546f231fbdSstefano_zampini   PetscErrorCode     ierr;
6556f231fbdSstefano_zampini   PetscFunctionBegin;
6566f231fbdSstefano_zampini   ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
6576f231fbdSstefano_zampini   ierr                   = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
6586f231fbdSstefano_zampini   (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_AIJ_wHYPRE;
659613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
660613e5ff0Sstefano_zampini }
661613e5ff0Sstefano_zampini 
6624cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_HYPRE(Mat A,Mat P,Mat C)
663613e5ff0Sstefano_zampini {
6644cc28894Sstefano_zampini   Mat                B;
6654cc28894Sstefano_zampini   Mat_HYPRE          *hP;
666613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *hA = NULL,*Pparcsr,*ptapparcsr;
667613e5ff0Sstefano_zampini   HYPRE_Int          type;
668613e5ff0Sstefano_zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
6694cc28894Sstefano_zampini   PetscBool          ishypre;
670613e5ff0Sstefano_zampini   PetscErrorCode     ierr;
671613e5ff0Sstefano_zampini 
672613e5ff0Sstefano_zampini   PetscFunctionBegin;
6734cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
6744cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(comm,PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
6754cc28894Sstefano_zampini   hP = (Mat_HYPRE*)P->data;
676613e5ff0Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
677613e5ff0Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
678613e5ff0Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
679613e5ff0Sstefano_zampini 
680613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
681613e5ff0Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,hA,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
682613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
683225daaf8SStefano Zampini 
6844cc28894Sstefano_zampini   /* create temporary matrix and merge to C */
6854cc28894Sstefano_zampini   ierr = MatCreateFromParCSR(ptapparcsr,((PetscObject)C)->type_name,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
6864cc28894Sstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
6874cc28894Sstefano_zampini   PetscFunctionReturn(0);
6884cc28894Sstefano_zampini }
6894cc28894Sstefano_zampini 
6904cc28894Sstefano_zampini static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
6914cc28894Sstefano_zampini {
6924cc28894Sstefano_zampini   PetscErrorCode ierr;
6934cc28894Sstefano_zampini 
6944cc28894Sstefano_zampini   PetscFunctionBegin;
6954cc28894Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
6964cc28894Sstefano_zampini     const char *deft = MATAIJ;
6974cc28894Sstefano_zampini     char       type[256];
6984cc28894Sstefano_zampini     PetscBool  flg;
6994cc28894Sstefano_zampini 
7004cc28894Sstefano_zampini     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
7014cc28894Sstefano_zampini     ierr = PetscOptionsFList("-matptap_hypre_outtype","Matrix type",NULL,MatList,deft,type,256,&flg);CHKERRQ(ierr);
7024cc28894Sstefano_zampini     ierr = PetscOptionsEnd();CHKERRQ(ierr);
7034cc28894Sstefano_zampini     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7044cc28894Sstefano_zampini     ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
7054cc28894Sstefano_zampini     if (flg) {
7064cc28894Sstefano_zampini       ierr = MatSetType(*C,type);CHKERRQ(ierr);
7074cc28894Sstefano_zampini     } else {
7084cc28894Sstefano_zampini       ierr = MatSetType(*C,deft);CHKERRQ(ierr);
7094cc28894Sstefano_zampini     }
7104cc28894Sstefano_zampini     (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_HYPRE;
7114cc28894Sstefano_zampini     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7124cc28894Sstefano_zampini   }
7134cc28894Sstefano_zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
7144cc28894Sstefano_zampini   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
715613e5ff0Sstefano_zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
716c1a070e6SStefano Zampini   PetscFunctionReturn(0);
717c1a070e6SStefano Zampini }
718c1a070e6SStefano Zampini 
7194cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_HYPRE_HYPRE(Mat A,Mat P,Mat C)
7204cc28894Sstefano_zampini {
7214cc28894Sstefano_zampini   Mat                B;
7224cc28894Sstefano_zampini   hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr;
7234cc28894Sstefano_zampini   Mat_HYPRE          *hA,*hP;
7244cc28894Sstefano_zampini   PetscBool          ishypre;
7254cc28894Sstefano_zampini   HYPRE_Int          type;
7264cc28894Sstefano_zampini   PetscErrorCode     ierr;
7274cc28894Sstefano_zampini 
7284cc28894Sstefano_zampini   PetscFunctionBegin;
7294cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
7304cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)P),PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
7314cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr);
7324cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE);
7334cc28894Sstefano_zampini   hA = (Mat_HYPRE*)A->data;
7344cc28894Sstefano_zampini   hP = (Mat_HYPRE*)P->data;
7354cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
7364cc28894Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
7374cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
7384cc28894Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
7394cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
7404cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
7414cc28894Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
7424cc28894Sstefano_zampini   ierr = MatCreateFromParCSR(ptapparcsr,MATHYPRE,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
7434cc28894Sstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
7444cc28894Sstefano_zampini   PetscFunctionReturn(0);
7454cc28894Sstefano_zampini }
7464cc28894Sstefano_zampini 
747cd8bc7baSStefano Zampini static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
748cd8bc7baSStefano Zampini {
749cd8bc7baSStefano Zampini   PetscErrorCode ierr;
750cd8bc7baSStefano Zampini 
751cd8bc7baSStefano Zampini   PetscFunctionBegin;
7524cc28894Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
7534cc28894Sstefano_zampini     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7544cc28894Sstefano_zampini     ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
7554cc28894Sstefano_zampini     ierr                   = MatSetType(*C,MATHYPRE);CHKERRQ(ierr);
7564cc28894Sstefano_zampini     (*C)->ops->ptapnumeric = MatPtAPNumeric_HYPRE_HYPRE;
7574cc28894Sstefano_zampini     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7584cc28894Sstefano_zampini   }
759613e5ff0Sstefano_zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
7604cc28894Sstefano_zampini   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
761613e5ff0Sstefano_zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
762cd8bc7baSStefano Zampini   PetscFunctionReturn(0);
763cd8bc7baSStefano Zampini }
764cd8bc7baSStefano Zampini 
765d501dc42Sstefano_zampini /* calls hypre_ParMatmul
766d501dc42Sstefano_zampini    hypre_ParMatMul uses hypre_ParMatrixCreate with the communicator of hA
7673dad0653Sstefano_zampini    hypre_ParMatrixCreate does not duplicate the communicator
7683dad0653Sstefano_zampini    It looks like we don't need to have the diagonal entries
7693dad0653Sstefano_zampini    ordered first in the rows of the diagonal part
7703dad0653Sstefano_zampini    for boomerAMGBuildCoarseOperator to work */
771d501dc42Sstefano_zampini static PetscErrorCode MatHYPRE_ParCSR_MatMatMult(hypre_ParCSRMatrix *hA, hypre_ParCSRMatrix *hB, hypre_ParCSRMatrix **hAB)
772d501dc42Sstefano_zampini {
773d501dc42Sstefano_zampini   PetscFunctionBegin;
774d501dc42Sstefano_zampini   PetscStackPush("hypre_ParMatmul");
775d501dc42Sstefano_zampini   *hAB = hypre_ParMatmul(hA,hB);
776d501dc42Sstefano_zampini   PetscStackPop;
777d501dc42Sstefano_zampini   PetscFunctionReturn(0);
778d501dc42Sstefano_zampini }
779d501dc42Sstefano_zampini 
7805e5acdf2Sstefano_zampini static PetscErrorCode MatMatMultNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C)
7815e5acdf2Sstefano_zampini {
7825e5acdf2Sstefano_zampini   Mat                D;
783d501dc42Sstefano_zampini   hypre_ParCSRMatrix *hA,*hB,*hAB = NULL;
7845e5acdf2Sstefano_zampini   PetscErrorCode     ierr;
7855e5acdf2Sstefano_zampini 
7865e5acdf2Sstefano_zampini   PetscFunctionBegin;
7875e5acdf2Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
7885e5acdf2Sstefano_zampini   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
789d501dc42Sstefano_zampini   ierr = MatHYPRE_ParCSR_MatMatMult(hA,hB,&hAB);CHKERRQ(ierr);
7905e5acdf2Sstefano_zampini   ierr = MatCreateFromParCSR(hAB,MATAIJ,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
7915e5acdf2Sstefano_zampini   ierr = MatHeaderMerge(C,&D);CHKERRQ(ierr);
7925e5acdf2Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
7935e5acdf2Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
7945e5acdf2Sstefano_zampini   PetscFunctionReturn(0);
7955e5acdf2Sstefano_zampini }
7965e5acdf2Sstefano_zampini 
7975e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat B,PetscReal fill,Mat *C)
7985e5acdf2Sstefano_zampini {
7995e5acdf2Sstefano_zampini   PetscErrorCode ierr;
80020e1dc0dSstefano_zampini 
8015e5acdf2Sstefano_zampini   PetscFunctionBegin;
8025e5acdf2Sstefano_zampini   ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
8035e5acdf2Sstefano_zampini   ierr                      = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
8045e5acdf2Sstefano_zampini   (*C)->ops->matmultnumeric = MatMatMultNumeric_AIJ_AIJ_wHYPRE;
8055e5acdf2Sstefano_zampini   PetscFunctionReturn(0);
8065e5acdf2Sstefano_zampini }
8075e5acdf2Sstefano_zampini 
808d501dc42Sstefano_zampini static PetscErrorCode MatMatMultNumeric_HYPRE_HYPRE(Mat A,Mat B,Mat C)
809d501dc42Sstefano_zampini {
810d501dc42Sstefano_zampini   Mat                D;
811d501dc42Sstefano_zampini   hypre_ParCSRMatrix *Aparcsr,*Bparcsr,*ABparcsr = NULL;
812d501dc42Sstefano_zampini   Mat_HYPRE          *hA,*hB;
813d501dc42Sstefano_zampini   PetscBool          ishypre;
814d501dc42Sstefano_zampini   HYPRE_Int          type;
815d501dc42Sstefano_zampini   PetscErrorCode     ierr;
816d501dc42Sstefano_zampini 
817d501dc42Sstefano_zampini   PetscFunctionBegin;
818d501dc42Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)B,MATHYPRE,&ishypre);CHKERRQ(ierr);
819d501dc42Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_USER,"B should be of type %s",MATHYPRE);
820d501dc42Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr);
821d501dc42Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE);
822d501dc42Sstefano_zampini   hA = (Mat_HYPRE*)A->data;
823d501dc42Sstefano_zampini   hB = (Mat_HYPRE*)B->data;
824d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
825d501dc42Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
826d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hB->ij,&type));
827d501dc42Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
828d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
829d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hB->ij,(void**)&Bparcsr));
830d501dc42Sstefano_zampini   ierr = MatHYPRE_ParCSR_MatMatMult(Aparcsr,Bparcsr,&ABparcsr);CHKERRQ(ierr);
831d501dc42Sstefano_zampini   ierr = MatCreateFromParCSR(ABparcsr,MATHYPRE,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
832d501dc42Sstefano_zampini   /* need to use HeaderReplace because HeaderMerge messes up with the communicator */
833d501dc42Sstefano_zampini   ierr = MatHeaderReplace(C,&D);CHKERRQ(ierr);
834d501dc42Sstefano_zampini   C->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE;
835d501dc42Sstefano_zampini   PetscFunctionReturn(0);
836d501dc42Sstefano_zampini }
837d501dc42Sstefano_zampini 
838d501dc42Sstefano_zampini static PetscErrorCode MatMatMult_HYPRE_HYPRE(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
839d501dc42Sstefano_zampini {
840d501dc42Sstefano_zampini   PetscErrorCode ierr;
841d501dc42Sstefano_zampini 
842d501dc42Sstefano_zampini   PetscFunctionBegin;
843d501dc42Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
844d501dc42Sstefano_zampini     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
845d501dc42Sstefano_zampini     ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
846d501dc42Sstefano_zampini     ierr                      = MatSetType(*C,MATHYPRE);CHKERRQ(ierr);
847d501dc42Sstefano_zampini     (*C)->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE;
848d501dc42Sstefano_zampini     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
849d501dc42Sstefano_zampini   }
850d501dc42Sstefano_zampini   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
851d501dc42Sstefano_zampini   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
852d501dc42Sstefano_zampini   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
853d501dc42Sstefano_zampini   PetscFunctionReturn(0);
854d501dc42Sstefano_zampini }
855d501dc42Sstefano_zampini 
8563dad0653Sstefano_zampini PETSC_INTERN PetscErrorCode MatTransposeMatMatMultNumeric_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,Mat D)
85720e1dc0dSstefano_zampini {
85820e1dc0dSstefano_zampini   Mat                E;
85920e1dc0dSstefano_zampini   hypre_ParCSRMatrix *hA,*hB,*hC,*hABC;
86020e1dc0dSstefano_zampini   PetscErrorCode     ierr;
86120e1dc0dSstefano_zampini 
86220e1dc0dSstefano_zampini   PetscFunctionBegin;
86320e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
86420e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
86520e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(C,&hC);CHKERRQ(ierr);
86620e1dc0dSstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(hA,hB,hC,&hABC);CHKERRQ(ierr);
86720e1dc0dSstefano_zampini   ierr = MatCreateFromParCSR(hABC,MATAIJ,PETSC_OWN_POINTER,&E);CHKERRQ(ierr);
86820e1dc0dSstefano_zampini   ierr = MatHeaderMerge(D,&E);CHKERRQ(ierr);
86920e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
87020e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
87120e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(C,&hC);CHKERRQ(ierr);
87220e1dc0dSstefano_zampini   PetscFunctionReturn(0);
87320e1dc0dSstefano_zampini }
87420e1dc0dSstefano_zampini 
8753dad0653Sstefano_zampini PETSC_INTERN PetscErrorCode MatTransposeMatMatMultSymbolic_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,PetscReal fill,Mat *D)
87620e1dc0dSstefano_zampini {
87720e1dc0dSstefano_zampini   PetscErrorCode ierr;
87820e1dc0dSstefano_zampini 
87920e1dc0dSstefano_zampini   PetscFunctionBegin;
88020e1dc0dSstefano_zampini   ierr = MatCreate(PetscObjectComm((PetscObject)A),D);CHKERRQ(ierr);
88120e1dc0dSstefano_zampini   ierr = MatSetType(*D,MATAIJ);CHKERRQ(ierr);
88220e1dc0dSstefano_zampini   PetscFunctionReturn(0);
88320e1dc0dSstefano_zampini }
88420e1dc0dSstefano_zampini 
885ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
88663c07aadSStefano Zampini {
88763c07aadSStefano Zampini   PetscErrorCode ierr;
88863c07aadSStefano Zampini 
88963c07aadSStefano Zampini   PetscFunctionBegin;
89063c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
89163c07aadSStefano Zampini   PetscFunctionReturn(0);
89263c07aadSStefano Zampini }
89363c07aadSStefano Zampini 
894ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
89563c07aadSStefano Zampini {
89663c07aadSStefano Zampini   PetscErrorCode ierr;
89763c07aadSStefano Zampini 
89863c07aadSStefano Zampini   PetscFunctionBegin;
89963c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
90063c07aadSStefano Zampini   PetscFunctionReturn(0);
90163c07aadSStefano Zampini }
90263c07aadSStefano Zampini 
903ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
90463c07aadSStefano Zampini {
90563c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
90663c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
90763c07aadSStefano Zampini   hypre_ParVector    *hx,*hy;
90863c07aadSStefano Zampini   PetscScalar        *ax,*ay,*sax,*say;
90963c07aadSStefano Zampini   PetscErrorCode     ierr;
91063c07aadSStefano Zampini 
91163c07aadSStefano Zampini   PetscFunctionBegin;
91263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
91363c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
91463c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
91563c07aadSStefano Zampini   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
91663c07aadSStefano Zampini   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
91763c07aadSStefano Zampini   if (trans) {
91858968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ay,say);
91958968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax);
92063c07aadSStefano Zampini     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
92158968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,say,ay);
92258968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax);
92363c07aadSStefano Zampini   } else {
92458968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax);
92558968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ay,say);
92663c07aadSStefano Zampini     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
92758968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax);
92858968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,say,ay);
92963c07aadSStefano Zampini   }
93063c07aadSStefano Zampini   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
93163c07aadSStefano Zampini   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
93263c07aadSStefano Zampini   PetscFunctionReturn(0);
93363c07aadSStefano Zampini }
93463c07aadSStefano Zampini 
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 
958ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A)
95963c07aadSStefano Zampini {
9604ec6421dSstefano_zampini   PetscErrorCode ierr;
9614ec6421dSstefano_zampini 
9624ec6421dSstefano_zampini   PetscFunctionBegin;
9634ec6421dSstefano_zampini   ierr = MatHYPRESetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr);
9644ec6421dSstefano_zampini   PetscFunctionReturn(0);
9654ec6421dSstefano_zampini }
9664ec6421dSstefano_zampini 
9674ec6421dSstefano_zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
9684ec6421dSstefano_zampini {
96963c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
97063c07aadSStefano Zampini   Vec                x,b;
97163c07aadSStefano Zampini   PetscErrorCode     ierr;
97263c07aadSStefano Zampini 
97363c07aadSStefano Zampini   PetscFunctionBegin;
9744ec6421dSstefano_zampini   if (mode == MAT_FLUSH_ASSEMBLY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_FLUSH_ASSEMBLY currently not supported with MATHYPRE");
9754ec6421dSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
97663c07aadSStefano Zampini   if (hA->x) PetscFunctionReturn(0);
97763c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
97863c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
97963c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
98063c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
98163c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
98263c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
98363c07aadSStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
98463c07aadSStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
98563c07aadSStefano Zampini   PetscFunctionReturn(0);
98663c07aadSStefano Zampini }
98763c07aadSStefano Zampini 
988d975228cSstefano_zampini #define MATHYPRE_SCRATCH 2048
989d975228cSstefano_zampini 
990d975228cSstefano_zampini PetscErrorCode MatSetValues_HYPRE(Mat A, PetscInt nr, const PetscInt rows[], PetscInt nc, const PetscInt cols[], const PetscScalar v[], InsertMode ins)
991d975228cSstefano_zampini {
992d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
993d975228cSstefano_zampini   PetscScalar        *vals = (PetscScalar *)v;
994d975228cSstefano_zampini   PetscScalar        sscr[MATHYPRE_SCRATCH];
9957d968826Sstefano_zampini   HYPRE_Int          cscr[2][MATHYPRE_SCRATCH];
9967d968826Sstefano_zampini   HYPRE_Int          i,nzc;
997d975228cSstefano_zampini   PetscErrorCode     ierr;
998d975228cSstefano_zampini 
999d975228cSstefano_zampini   PetscFunctionBegin;
1000d975228cSstefano_zampini   for (i=0,nzc=0;i<nc;i++) {
1001d975228cSstefano_zampini     if (cols[i] >= 0) {
1002d975228cSstefano_zampini       cscr[0][nzc  ] = cols[i];
1003d975228cSstefano_zampini       cscr[1][nzc++] = i;
1004d975228cSstefano_zampini     }
1005d975228cSstefano_zampini   }
1006d975228cSstefano_zampini   if (!nzc) PetscFunctionReturn(0);
1007d975228cSstefano_zampini 
1008d975228cSstefano_zampini   if (ins == ADD_VALUES) {
1009d975228cSstefano_zampini     for (i=0;i<nr;i++) {
1010e3977e59Sstefano_zampini       if (rows[i] >= 0 && nzc) {
1011d975228cSstefano_zampini         PetscInt j;
1012d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
10137d968826Sstefano_zampini         PetscStackCallStandard(HYPRE_IJMatrixAddToValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
1014d975228cSstefano_zampini       }
1015d975228cSstefano_zampini       vals += nc;
1016d975228cSstefano_zampini     }
1017d975228cSstefano_zampini   } else { /* INSERT_VALUES */
1018d975228cSstefano_zampini #if defined(PETSC_USE_DEBUG)
1019d975228cSstefano_zampini     /* Insert values cannot be used to insert offproc entries */
1020d975228cSstefano_zampini     PetscInt rst,ren;
1021d975228cSstefano_zampini     ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
1022d975228cSstefano_zampini     for (i=0;i<nr;i++)
1023d975228cSstefano_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");
1024d975228cSstefano_zampini #endif
1025d975228cSstefano_zampini     for (i=0;i<nr;i++) {
1026e3977e59Sstefano_zampini       if (rows[i] >= 0 && nzc) {
1027d975228cSstefano_zampini         PetscInt j;
1028d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
10297d968826Sstefano_zampini         PetscStackCallStandard(HYPRE_IJMatrixSetValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
1030d975228cSstefano_zampini       }
1031d975228cSstefano_zampini       vals += nc;
1032d975228cSstefano_zampini     }
1033d975228cSstefano_zampini   }
1034d975228cSstefano_zampini   PetscFunctionReturn(0);
1035d975228cSstefano_zampini }
1036d975228cSstefano_zampini 
1037d975228cSstefano_zampini static PetscErrorCode MatHYPRESetPreallocation_HYPRE(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1038d975228cSstefano_zampini {
1039d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
10407d968826Sstefano_zampini   HYPRE_Int          *hdnnz,*honnz;
104106a29025Sstefano_zampini   PetscInt           i,rs,re,cs,ce,bs;
1042d975228cSstefano_zampini   PetscMPIInt        size;
1043d975228cSstefano_zampini   PetscErrorCode     ierr;
1044d975228cSstefano_zampini 
1045d975228cSstefano_zampini   PetscFunctionBegin;
104606a29025Sstefano_zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1047d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1048d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1049d975228cSstefano_zampini   rs   = A->rmap->rstart;
1050d975228cSstefano_zampini   re   = A->rmap->rend;
1051d975228cSstefano_zampini   cs   = A->cmap->rstart;
1052d975228cSstefano_zampini   ce   = A->cmap->rend;
1053d975228cSstefano_zampini   if (!hA->ij) {
1054d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rs,re-1,cs,ce-1,&hA->ij));
1055d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1056d975228cSstefano_zampini   } else {
10577d968826Sstefano_zampini     HYPRE_Int hrs,hre,hcs,hce;
1058d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixGetLocalRange,(hA->ij,&hrs,&hre,&hcs,&hce));
1059d975228cSstefano_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);
1060d975228cSstefano_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);
1061d975228cSstefano_zampini   }
1062d975228cSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1063d975228cSstefano_zampini 
106406a29025Sstefano_zampini   if (dnz == PETSC_DEFAULT || dnz == PETSC_DECIDE) dnz = 10*bs;
106506a29025Sstefano_zampini   if (onz == PETSC_DEFAULT || onz == PETSC_DECIDE) onz = 10*bs;
106606a29025Sstefano_zampini 
1067d975228cSstefano_zampini   if (!dnnz) {
1068d975228cSstefano_zampini     ierr = PetscMalloc1(A->rmap->n,&hdnnz);CHKERRQ(ierr);
1069d975228cSstefano_zampini     for (i=0;i<A->rmap->n;i++) hdnnz[i] = dnz;
1070d975228cSstefano_zampini   } else {
10717d968826Sstefano_zampini     hdnnz = (HYPRE_Int*)dnnz;
1072d975228cSstefano_zampini   }
1073d975228cSstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1074d975228cSstefano_zampini   if (size > 1) {
1075d975228cSstefano_zampini     if (!onnz) {
1076d975228cSstefano_zampini       ierr = PetscMalloc1(A->rmap->n,&honnz);CHKERRQ(ierr);
1077d975228cSstefano_zampini       for (i=0;i<A->rmap->n;i++) honnz[i] = onz;
1078d975228cSstefano_zampini     } else {
10797d968826Sstefano_zampini       honnz = (HYPRE_Int*)onnz;
1080d975228cSstefano_zampini     }
1081d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(hA->ij,hdnnz,honnz));
1082d975228cSstefano_zampini   } else {
1083d975228cSstefano_zampini     honnz = NULL;
1084d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetRowSizes,(hA->ij,hdnnz));
1085d975228cSstefano_zampini   }
1086d975228cSstefano_zampini   if (!dnnz) {
1087d975228cSstefano_zampini     ierr = PetscFree(hdnnz);CHKERRQ(ierr);
1088d975228cSstefano_zampini   }
1089d975228cSstefano_zampini   if (!onnz && honnz) {
1090d975228cSstefano_zampini     ierr = PetscFree(honnz);CHKERRQ(ierr);
1091d975228cSstefano_zampini   }
109206a29025Sstefano_zampini   A->preallocated = PETSC_TRUE;
1093d975228cSstefano_zampini 
1094d975228cSstefano_zampini   /* SetDiagOffdSizes sets hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0 */
1095d975228cSstefano_zampini   {
1096d975228cSstefano_zampini     hypre_AuxParCSRMatrix *aux_matrix;
1097d975228cSstefano_zampini     aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij);
1098d975228cSstefano_zampini     hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1;
1099d975228cSstefano_zampini   }
1100d975228cSstefano_zampini   PetscFunctionReturn(0);
1101d975228cSstefano_zampini }
1102d975228cSstefano_zampini 
1103d975228cSstefano_zampini /*@C
1104d975228cSstefano_zampini    MatHYPRESetPreallocation - Preallocates memory for a sparse parallel matrix in HYPRE IJ format
1105d975228cSstefano_zampini 
1106d975228cSstefano_zampini    Collective on Mat
1107d975228cSstefano_zampini 
1108d975228cSstefano_zampini    Input Parameters:
1109d975228cSstefano_zampini +  A - the matrix
1110d975228cSstefano_zampini .  dnz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1111d975228cSstefano_zampini           (same value is used for all local rows)
1112d975228cSstefano_zampini .  dnnz - array containing the number of nonzeros in the various rows of the
1113d975228cSstefano_zampini           DIAGONAL portion of the local submatrix (possibly different for each row)
1114d975228cSstefano_zampini           or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure.
1115d975228cSstefano_zampini           The size of this array is equal to the number of local rows, i.e 'm'.
1116d975228cSstefano_zampini           For matrices that will be factored, you must leave room for (and set)
1117d975228cSstefano_zampini           the diagonal entry even if it is zero.
1118d975228cSstefano_zampini .  onz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1119d975228cSstefano_zampini           submatrix (same value is used for all local rows).
1120d975228cSstefano_zampini -  onnz - array containing the number of nonzeros in the various rows of the
1121d975228cSstefano_zampini           OFF-DIAGONAL portion of the local submatrix (possibly different for
1122d975228cSstefano_zampini           each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero
1123d975228cSstefano_zampini           structure. The size of this array is equal to the number
1124d975228cSstefano_zampini           of local rows, i.e 'm'.
1125d975228cSstefano_zampini 
1126d975228cSstefano_zampini    Notes: If the *nnz parameter is given then the *nz parameter is ignored; for sequential matrices, onz and onnz are ignored.
1127d975228cSstefano_zampini 
1128d975228cSstefano_zampini    Level: intermediate
1129d975228cSstefano_zampini 
1130d975228cSstefano_zampini .keywords: matrix, aij, compressed row, sparse, parallel
1131d975228cSstefano_zampini 
1132d975228cSstefano_zampini .seealso: MatCreate(), MatMPIAIJSetPreallocation, MATHYPRE
1133d975228cSstefano_zampini @*/
1134d975228cSstefano_zampini PetscErrorCode MatHYPRESetPreallocation(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1135d975228cSstefano_zampini {
1136d975228cSstefano_zampini   PetscErrorCode ierr;
1137d975228cSstefano_zampini 
1138d975228cSstefano_zampini   PetscFunctionBegin;
1139d975228cSstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1140d975228cSstefano_zampini   PetscValidType(A,1);
1141d975228cSstefano_zampini   ierr = PetscTryMethod(A,"MatHYPRESetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(A,dnz,dnnz,onz,onnz));CHKERRQ(ierr);
1142d975228cSstefano_zampini   PetscFunctionReturn(0);
1143d975228cSstefano_zampini }
1144d975228cSstefano_zampini 
1145225daaf8SStefano Zampini /*
1146225daaf8SStefano Zampini    MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix
1147225daaf8SStefano Zampini 
1148225daaf8SStefano Zampini    Collective
1149225daaf8SStefano Zampini 
1150225daaf8SStefano Zampini    Input Parameters:
1151225daaf8SStefano Zampini +  vparcsr  - the pointer to the hypre_ParCSRMatrix
1152bb4689ddSStefano Zampini .  mtype    - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported.
1153225daaf8SStefano Zampini -  copymode - PETSc copying options
1154225daaf8SStefano Zampini 
1155225daaf8SStefano Zampini    Output Parameter:
1156225daaf8SStefano Zampini .  A  - the matrix
1157225daaf8SStefano Zampini 
1158225daaf8SStefano Zampini    Level: intermediate
1159225daaf8SStefano Zampini 
1160225daaf8SStefano Zampini .seealso: MatHYPRE, PetscCopyMode
1161225daaf8SStefano Zampini */
1162dd9c0a25Sstefano_zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(hypre_ParCSRMatrix *vparcsr, MatType mtype, PetscCopyMode copymode, Mat* A)
1163978814f1SStefano Zampini {
1164225daaf8SStefano Zampini   Mat                   T;
1165978814f1SStefano Zampini   Mat_HYPRE             *hA;
116663c07aadSStefano Zampini   hypre_ParCSRMatrix    *parcsr;
1167978814f1SStefano Zampini   MPI_Comm              comm;
1168978814f1SStefano Zampini   PetscInt              rstart,rend,cstart,cend,M,N;
1169bb4689ddSStefano Zampini   PetscBool             isseqaij,ismpiaij,isaij,ishyp,isis;
1170978814f1SStefano Zampini   PetscErrorCode        ierr;
1171978814f1SStefano Zampini 
1172978814f1SStefano Zampini   PetscFunctionBegin;
1173978814f1SStefano Zampini   parcsr = (hypre_ParCSRMatrix *)vparcsr;
1174978814f1SStefano Zampini   comm   = hypre_ParCSRMatrixComm(parcsr);
1175225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
1176225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
1177225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr);
1178225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr);
11798cfe8d00SStefano Zampini   ierr   = PetscStrcmp(mtype,MATIS,&isis);CHKERRQ(ierr);
1180225daaf8SStefano Zampini   isaij  = (PetscBool)(isseqaij || ismpiaij || isaij);
1181bb4689ddSStefano 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);
1182225daaf8SStefano Zampini   if (ishyp && copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES");
1183978814f1SStefano Zampini 
1184978814f1SStefano Zampini   /* access ParCSRMatrix */
1185978814f1SStefano Zampini   rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr);
1186978814f1SStefano Zampini   rend   = hypre_ParCSRMatrixLastRowIndex(parcsr);
1187978814f1SStefano Zampini   cstart = hypre_ParCSRMatrixFirstColDiag(parcsr);
1188978814f1SStefano Zampini   cend   = hypre_ParCSRMatrixLastColDiag(parcsr);
1189978814f1SStefano Zampini   M      = hypre_ParCSRMatrixGlobalNumRows(parcsr);
1190978814f1SStefano Zampini   N      = hypre_ParCSRMatrixGlobalNumCols(parcsr);
1191978814f1SStefano Zampini 
1192fa92c42cSstefano_zampini   /* fix for empty local rows/columns */
1193fa92c42cSstefano_zampini   if (rend < rstart) rend = rstart;
1194fa92c42cSstefano_zampini   if (cend < cstart) cend = cstart;
1195fa92c42cSstefano_zampini 
1196978814f1SStefano Zampini   /* create PETSc matrix with MatHYPRE */
1197225daaf8SStefano Zampini   ierr = MatCreate(comm,&T);CHKERRQ(ierr);
1198225daaf8SStefano Zampini   ierr = MatSetSizes(T,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr);
1199225daaf8SStefano Zampini   ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr);
1200225daaf8SStefano Zampini   hA   = (Mat_HYPRE*)(T->data);
1201978814f1SStefano Zampini 
1202978814f1SStefano Zampini   /* create HYPRE_IJMatrix */
1203978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
1204978814f1SStefano Zampini 
1205978814f1SStefano Zampini   /* set ParCSR object */
1206978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1207978814f1SStefano Zampini   hypre_IJMatrixObject(hA->ij) = parcsr;
12084ec6421dSstefano_zampini   T->preallocated = PETSC_TRUE;
1209978814f1SStefano Zampini 
1210978814f1SStefano Zampini   /* set assembled flag */
1211978814f1SStefano Zampini   hypre_IJMatrixAssembleFlag(hA->ij) = 1;
1212978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1213225daaf8SStefano Zampini   if (ishyp) {
12146d2a658fSstefano_zampini     PetscMPIInt myid = 0;
12156d2a658fSstefano_zampini 
12166d2a658fSstefano_zampini     /* make sure we always have row_starts and col_starts available */
12176d2a658fSstefano_zampini     if (HYPRE_AssumedPartitionCheck()) {
12186d2a658fSstefano_zampini       ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
12196d2a658fSstefano_zampini     }
12206d2a658fSstefano_zampini     if (!hypre_ParCSRMatrixOwnsColStarts(parcsr)) {
12216d2a658fSstefano_zampini       PetscLayout map;
12226d2a658fSstefano_zampini 
12236d2a658fSstefano_zampini       ierr = MatGetLayouts(T,NULL,&map);CHKERRQ(ierr);
12246d2a658fSstefano_zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
12251c92d2d0Sstefano_zampini       hypre_ParCSRMatrixColStarts(parcsr) = (HYPRE_Int*)(map->range + myid);
12266d2a658fSstefano_zampini     }
12276d2a658fSstefano_zampini     if (!hypre_ParCSRMatrixOwnsRowStarts(parcsr)) {
12286d2a658fSstefano_zampini       PetscLayout map;
12296d2a658fSstefano_zampini 
12306d2a658fSstefano_zampini       ierr = MatGetLayouts(T,&map,NULL);CHKERRQ(ierr);
12316d2a658fSstefano_zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
12321c92d2d0Sstefano_zampini       hypre_ParCSRMatrixRowStarts(parcsr) = (HYPRE_Int*)(map->range + myid);
12336d2a658fSstefano_zampini     }
1234978814f1SStefano Zampini     /* prevent from freeing the pointer */
1235978814f1SStefano Zampini     if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE;
1236225daaf8SStefano Zampini     *A   = T;
1237978814f1SStefano Zampini     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1238978814f1SStefano Zampini     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1239bb4689ddSStefano Zampini   } else if (isaij) {
1240bb4689ddSStefano Zampini     if (copymode != PETSC_OWN_POINTER) {
1241225daaf8SStefano Zampini       /* prevent from freeing the pointer */
1242225daaf8SStefano Zampini       hA->inner_free = PETSC_FALSE;
1243225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
1244225daaf8SStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
1245225daaf8SStefano Zampini     } else { /* AIJ return type with PETSC_OWN_POINTER */
1246225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1247225daaf8SStefano Zampini       *A   = T;
1248225daaf8SStefano Zampini     }
1249bb4689ddSStefano Zampini   } else if (isis) {
12508cfe8d00SStefano Zampini     ierr = MatConvert_HYPRE_IS(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
12518cfe8d00SStefano Zampini     if (copymode != PETSC_OWN_POINTER) hA->inner_free = PETSC_FALSE;
1252bb4689ddSStefano Zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
1253bb4689ddSStefano Zampini   }
1254978814f1SStefano Zampini   PetscFunctionReturn(0);
1255978814f1SStefano Zampini }
1256978814f1SStefano Zampini 
1257dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR_HYPRE(Mat A, hypre_ParCSRMatrix **parcsr)
1258dd9c0a25Sstefano_zampini {
1259dd9c0a25Sstefano_zampini   Mat_HYPRE*            hA = (Mat_HYPRE*)A->data;
1260dd9c0a25Sstefano_zampini   HYPRE_Int             type;
1261dd9c0a25Sstefano_zampini   PetscErrorCode        ierr;
1262dd9c0a25Sstefano_zampini 
1263dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1264dd9c0a25Sstefano_zampini   if (!hA->ij) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"HYPRE_IJMatrix not present");
1265dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
1266dd9c0a25Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"HYPRE_IJMatrix is not of type HYPRE_PARCSR");
1267dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)parcsr));
1268dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1269dd9c0a25Sstefano_zampini }
1270dd9c0a25Sstefano_zampini 
1271dd9c0a25Sstefano_zampini /*
1272dd9c0a25Sstefano_zampini    MatHYPREGetParCSR - Gets the pointer to the ParCSR matrix
1273dd9c0a25Sstefano_zampini 
1274dd9c0a25Sstefano_zampini    Not collective
1275dd9c0a25Sstefano_zampini 
1276dd9c0a25Sstefano_zampini    Input Parameters:
1277dd9c0a25Sstefano_zampini +  A  - the MATHYPRE object
1278dd9c0a25Sstefano_zampini 
1279dd9c0a25Sstefano_zampini    Output Parameter:
1280dd9c0a25Sstefano_zampini .  parcsr  - the pointer to the hypre_ParCSRMatrix
1281dd9c0a25Sstefano_zampini 
1282dd9c0a25Sstefano_zampini    Level: intermediate
1283dd9c0a25Sstefano_zampini 
1284dd9c0a25Sstefano_zampini .seealso: MatHYPRE, PetscCopyMode
1285dd9c0a25Sstefano_zampini */
1286dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR(Mat A, hypre_ParCSRMatrix **parcsr)
1287dd9c0a25Sstefano_zampini {
1288dd9c0a25Sstefano_zampini   PetscErrorCode ierr;
1289dd9c0a25Sstefano_zampini 
1290dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1291dd9c0a25Sstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1292dd9c0a25Sstefano_zampini   PetscValidType(A,1);
1293dd9c0a25Sstefano_zampini   ierr = PetscUseMethod(A,"MatHYPREGetParCSR_C",(Mat,hypre_ParCSRMatrix**),(A,parcsr));CHKERRQ(ierr);
1294dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1295dd9c0a25Sstefano_zampini }
1296dd9c0a25Sstefano_zampini 
1297*a055b5aaSBarry Smith /*MC
1298*a055b5aaSBarry Smith    MATHYPRE - MATHYPRE = "hypre" - A matrix type to be used for sequential and parallel sparse matrices
1299*a055b5aaSBarry Smith           based on the hypre IJ interface.
1300*a055b5aaSBarry Smith 
1301*a055b5aaSBarry Smith    Level: intermediate
1302*a055b5aaSBarry Smith 
1303*a055b5aaSBarry Smith .seealso: MatCreate()
1304*a055b5aaSBarry Smith M*/
1305*a055b5aaSBarry Smith 
130663c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
130763c07aadSStefano Zampini {
130863c07aadSStefano Zampini   Mat_HYPRE      *hB;
130963c07aadSStefano Zampini   PetscErrorCode ierr;
131063c07aadSStefano Zampini 
131163c07aadSStefano Zampini   PetscFunctionBegin;
131263c07aadSStefano Zampini   ierr           = PetscNewLog(B,&hB);CHKERRQ(ierr);
1313978814f1SStefano Zampini   hB->inner_free = PETSC_TRUE;
1314978814f1SStefano Zampini 
131563c07aadSStefano Zampini   B->data       = (void*)hB;
131663c07aadSStefano Zampini   B->rmap->bs   = 1;
131763c07aadSStefano Zampini   B->assembled  = PETSC_FALSE;
131863c07aadSStefano Zampini 
131963c07aadSStefano Zampini   B->ops->mult          = MatMult_HYPRE;
132063c07aadSStefano Zampini   B->ops->multtranspose = MatMultTranspose_HYPRE;
132163c07aadSStefano Zampini   B->ops->setup         = MatSetUp_HYPRE;
132263c07aadSStefano Zampini   B->ops->destroy       = MatDestroy_HYPRE;
132363c07aadSStefano Zampini   B->ops->assemblyend   = MatAssemblyEnd_HYPRE;
1324cd8bc7baSStefano Zampini   B->ops->ptap          = MatPtAP_HYPRE_HYPRE;
1325d501dc42Sstefano_zampini   B->ops->matmult       = MatMatMult_HYPRE_HYPRE;
1326d975228cSstefano_zampini   B->ops->setvalues     = MatSetValues_HYPRE;
132763c07aadSStefano Zampini 
132863c07aadSStefano Zampini   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
132963c07aadSStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
133063c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
13312df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr);
1332c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1333c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1334d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPRESetPreallocation_C",MatHYPRESetPreallocation_HYPRE);CHKERRQ(ierr);
1335dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPREGetParCSR_C",MatHYPREGetParCSR_HYPRE);CHKERRQ(ierr);
133663c07aadSStefano Zampini   PetscFunctionReturn(0);
133763c07aadSStefano Zampini }
133863c07aadSStefano Zampini 
1339225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void *ptr)
1340225daaf8SStefano Zampini {
1341225daaf8SStefano Zampini    PetscFunctionBegin;
1342225daaf8SStefano Zampini    hypre_TFree(ptr);
1343225daaf8SStefano Zampini    PetscFunctionReturn(0);
1344225daaf8SStefano Zampini }
1345