xref: /petsc/src/mat/impls/hypre/mhypre.c (revision 08defe43e687255440fe9c8a27f49f6d411c3f61)
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>
1468ec7858SStefano Zampini #include <_hypre_sstruct_ls.h>
1563c07aadSStefano Zampini 
1663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*);
1763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix);
1863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
1963c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
2063c07aadSStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool);
21225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void*);
22c69f721fSFande Kong PetscErrorCode MatSetValues_HYPRE(Mat, PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode ins);
2363c07aadSStefano Zampini 
2463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij)
2563c07aadSStefano Zampini {
2663c07aadSStefano Zampini   PetscErrorCode ierr;
2763c07aadSStefano Zampini   PetscInt       i,n_d,n_o;
2863c07aadSStefano Zampini   const PetscInt *ia_d,*ia_o;
2963c07aadSStefano Zampini   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
3063c07aadSStefano Zampini   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
3163c07aadSStefano Zampini 
3263c07aadSStefano Zampini   PetscFunctionBegin;
3363c07aadSStefano Zampini   if (A_d) { /* determine number of nonzero entries in local diagonal part */
3463c07aadSStefano Zampini     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
3563c07aadSStefano Zampini     if (done_d) {
3663c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr);
3763c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
3863c07aadSStefano Zampini         nnz_d[i] = ia_d[i+1] - ia_d[i];
3963c07aadSStefano Zampini       }
4063c07aadSStefano Zampini     }
4163c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
4263c07aadSStefano Zampini   }
4363c07aadSStefano Zampini   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
4463c07aadSStefano Zampini     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
4563c07aadSStefano Zampini     if (done_o) {
4663c07aadSStefano Zampini       ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr);
4763c07aadSStefano Zampini       for (i=0; i<n_o; i++) {
4863c07aadSStefano Zampini         nnz_o[i] = ia_o[i+1] - ia_o[i];
4963c07aadSStefano Zampini       }
5063c07aadSStefano Zampini     }
5163c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
5263c07aadSStefano Zampini   }
5363c07aadSStefano Zampini   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
5463c07aadSStefano Zampini     if (!done_o) { /* only diagonal part */
5563c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
5663c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
5763c07aadSStefano Zampini         nnz_o[i] = 0;
5863c07aadSStefano Zampini       }
5963c07aadSStefano Zampini     }
6063c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o));
6163c07aadSStefano Zampini     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
6263c07aadSStefano Zampini     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
6363c07aadSStefano Zampini   }
6463c07aadSStefano Zampini   PetscFunctionReturn(0);
6563c07aadSStefano Zampini }
6663c07aadSStefano Zampini 
6763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA)
6863c07aadSStefano Zampini {
6963c07aadSStefano Zampini   PetscErrorCode ierr;
7063c07aadSStefano Zampini   PetscInt       rstart,rend,cstart,cend;
7163c07aadSStefano Zampini 
7263c07aadSStefano Zampini   PetscFunctionBegin;
73d975228cSstefano_zampini   ierr   = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
74d975228cSstefano_zampini   ierr   = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
7563c07aadSStefano Zampini   rstart = A->rmap->rstart;
7663c07aadSStefano Zampini   rend   = A->rmap->rend;
7763c07aadSStefano Zampini   cstart = A->cmap->rstart;
7863c07aadSStefano Zampini   cend   = A->cmap->rend;
7963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
8063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
8163c07aadSStefano Zampini   {
8263c07aadSStefano Zampini     PetscBool      same;
8363c07aadSStefano Zampini     Mat            A_d,A_o;
8463c07aadSStefano Zampini     const PetscInt *colmap;
8563c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
8663c07aadSStefano Zampini     if (same) {
8763c07aadSStefano Zampini       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
8863c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
8963c07aadSStefano Zampini       PetscFunctionReturn(0);
9063c07aadSStefano Zampini     }
9163c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
9263c07aadSStefano Zampini     if (same) {
9363c07aadSStefano Zampini       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
9463c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
9563c07aadSStefano Zampini       PetscFunctionReturn(0);
9663c07aadSStefano Zampini     }
9763c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
9863c07aadSStefano Zampini     if (same) {
9963c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
10063c07aadSStefano Zampini       PetscFunctionReturn(0);
10163c07aadSStefano Zampini     }
10263c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
10363c07aadSStefano Zampini     if (same) {
10463c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
10563c07aadSStefano Zampini       PetscFunctionReturn(0);
10663c07aadSStefano Zampini     }
10763c07aadSStefano Zampini   }
10863c07aadSStefano Zampini   PetscFunctionReturn(0);
10963c07aadSStefano Zampini }
11063c07aadSStefano Zampini 
11163c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij)
11263c07aadSStefano Zampini {
11363c07aadSStefano Zampini   PetscErrorCode    ierr;
11463c07aadSStefano Zampini   PetscInt          i,rstart,rend,ncols,nr,nc;
11563c07aadSStefano Zampini   const PetscScalar *values;
11663c07aadSStefano Zampini   const PetscInt    *cols;
11763c07aadSStefano Zampini   PetscBool         flg;
11863c07aadSStefano Zampini 
11963c07aadSStefano Zampini   PetscFunctionBegin;
12063c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
12163c07aadSStefano Zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
12263c07aadSStefano Zampini   if (flg && nr == nc) {
12363c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
12463c07aadSStefano Zampini     PetscFunctionReturn(0);
12563c07aadSStefano Zampini   }
12663c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
12763c07aadSStefano Zampini   if (flg) {
12863c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
12963c07aadSStefano Zampini     PetscFunctionReturn(0);
13063c07aadSStefano Zampini   }
13163c07aadSStefano Zampini 
13263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
13363c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
13463c07aadSStefano Zampini   for (i=rstart; i<rend; i++) {
13563c07aadSStefano Zampini     ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
136e3977e59Sstefano_zampini     if (ncols) {
13763c07aadSStefano Zampini       PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
138e3977e59Sstefano_zampini     }
13963c07aadSStefano Zampini     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
14063c07aadSStefano Zampini   }
14163c07aadSStefano Zampini   PetscFunctionReturn(0);
14263c07aadSStefano Zampini }
14363c07aadSStefano Zampini 
14463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij)
14563c07aadSStefano Zampini {
14663c07aadSStefano Zampini   PetscErrorCode        ierr;
14763c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag = (Mat_SeqAIJ*)A->data;
14858968eb6SStefano Zampini   HYPRE_Int             type;
14963c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
15063c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
15163c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag;
15263c07aadSStefano Zampini 
15363c07aadSStefano Zampini   PetscFunctionBegin;
15463c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
155ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
156ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
157ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
15863c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
15963c07aadSStefano Zampini   /*
16063c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
16163c07aadSStefano Zampini   */
16263c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));
16363c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));
16463c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
165ea9daf28SStefano Zampini 
166ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
16763c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
16863c07aadSStefano Zampini   PetscFunctionReturn(0);
16963c07aadSStefano Zampini }
17063c07aadSStefano Zampini 
17163c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
17263c07aadSStefano Zampini {
17363c07aadSStefano Zampini   PetscErrorCode        ierr;
17463c07aadSStefano Zampini   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
17563c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag,*poffd;
17663c07aadSStefano Zampini   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
17758968eb6SStefano Zampini   HYPRE_Int             type;
17863c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
17963c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
18063c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag,*hoffd;
18163c07aadSStefano Zampini 
18263c07aadSStefano Zampini   PetscFunctionBegin;
18363c07aadSStefano Zampini   pdiag = (Mat_SeqAIJ*) pA->A->data;
18463c07aadSStefano Zampini   poffd = (Mat_SeqAIJ*) pA->B->data;
18563c07aadSStefano Zampini   /* cstart is only valid for square MPIAIJ layed out in the usual way */
18663c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
18763c07aadSStefano Zampini 
18863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
189ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
190ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
191ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
19263c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
19363c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
19463c07aadSStefano Zampini 
19563c07aadSStefano Zampini   /*
19663c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
19763c07aadSStefano Zampini   */
19863c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
19963c07aadSStefano Zampini   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
20063c07aadSStefano Zampini   jj  = (PetscInt*)hdiag->j;
20163c07aadSStefano Zampini   pjj = (PetscInt*)pdiag->j;
20263c07aadSStefano Zampini   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
20363c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
20463c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
20563c07aadSStefano Zampini   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
20663c07aadSStefano Zampini      If we hacked a hypre a bit more we might be able to avoid this step */
20763c07aadSStefano Zampini   jj  = (PetscInt*) hoffd->j;
20863c07aadSStefano Zampini   pjj = (PetscInt*) poffd->j;
20963c07aadSStefano Zampini   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
21063c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
21163c07aadSStefano Zampini 
212ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
21363c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
21463c07aadSStefano Zampini   PetscFunctionReturn(0);
21563c07aadSStefano Zampini }
21663c07aadSStefano Zampini 
2172df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B)
2182df22349SStefano Zampini {
2192df22349SStefano Zampini   Mat_HYPRE*             mhA = (Mat_HYPRE*)(A->data);
2202df22349SStefano Zampini   Mat                    lA;
2212df22349SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
2222df22349SStefano Zampini   IS                     is;
2232df22349SStefano Zampini   hypre_ParCSRMatrix     *hA;
2242df22349SStefano Zampini   hypre_CSRMatrix        *hdiag,*hoffd;
2252df22349SStefano Zampini   MPI_Comm               comm;
2262df22349SStefano Zampini   PetscScalar            *hdd,*hod,*aa,*data;
2277d968826Sstefano_zampini   HYPRE_Int              *col_map_offd,*hdi,*hdj,*hoi,*hoj;
2282df22349SStefano Zampini   PetscInt               *ii,*jj,*iptr,*jptr;
2292df22349SStefano Zampini   PetscInt               cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N;
23058968eb6SStefano Zampini   HYPRE_Int              type;
2312df22349SStefano Zampini   PetscErrorCode         ierr;
2322df22349SStefano Zampini 
2332df22349SStefano Zampini   PetscFunctionBegin;
234a1787963SStefano Zampini   comm = PetscObjectComm((PetscObject)A);
2352df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type));
2362df22349SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
2372df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA));
2382df22349SStefano Zampini   M     = hypre_ParCSRMatrixGlobalNumRows(hA);
2392df22349SStefano Zampini   N     = hypre_ParCSRMatrixGlobalNumCols(hA);
2402df22349SStefano Zampini   str   = hypre_ParCSRMatrixFirstRowIndex(hA);
2412df22349SStefano Zampini   stc   = hypre_ParCSRMatrixFirstColDiag(hA);
2422df22349SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(hA);
2432df22349SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(hA);
2442df22349SStefano Zampini   dr    = hypre_CSRMatrixNumRows(hdiag);
2452df22349SStefano Zampini   dc    = hypre_CSRMatrixNumCols(hdiag);
2462df22349SStefano Zampini   nnz   = hypre_CSRMatrixNumNonzeros(hdiag);
2472df22349SStefano Zampini   hdi   = hypre_CSRMatrixI(hdiag);
2482df22349SStefano Zampini   hdj   = hypre_CSRMatrixJ(hdiag);
2492df22349SStefano Zampini   hdd   = hypre_CSRMatrixData(hdiag);
2502df22349SStefano Zampini   oc    = hypre_CSRMatrixNumCols(hoffd);
2512df22349SStefano Zampini   nnz  += hypre_CSRMatrixNumNonzeros(hoffd);
2522df22349SStefano Zampini   hoi   = hypre_CSRMatrixI(hoffd);
2532df22349SStefano Zampini   hoj   = hypre_CSRMatrixJ(hoffd);
2542df22349SStefano Zampini   hod   = hypre_CSRMatrixData(hoffd);
2552df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
2562df22349SStefano Zampini     PetscInt *aux;
2572df22349SStefano Zampini 
2582df22349SStefano Zampini     /* generate l2g maps for rows and cols */
2592df22349SStefano Zampini     ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr);
2602df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
2612df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2622df22349SStefano Zampini     col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
2632df22349SStefano Zampini     ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr);
2642df22349SStefano Zampini     for (i=0; i<dc; i++) aux[i] = i+stc;
2652df22349SStefano Zampini     for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i];
2662df22349SStefano Zampini     ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
2672df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
2682df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2692df22349SStefano Zampini     /* create MATIS object */
2702df22349SStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
2712df22349SStefano Zampini     ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr);
2722df22349SStefano Zampini     ierr = MatSetType(*B,MATIS);CHKERRQ(ierr);
2732df22349SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr);
2742df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
2752df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
2762df22349SStefano Zampini 
2772df22349SStefano Zampini     /* allocate CSR for local matrix */
2782df22349SStefano Zampini     ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr);
2792df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr);
2802df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr);
2812df22349SStefano Zampini   } else {
2822df22349SStefano Zampini     PetscInt  nr;
2832df22349SStefano Zampini     PetscBool done;
2842df22349SStefano Zampini     ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr);
2852df22349SStefano Zampini     ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr);
2862df22349SStefano 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);
2872df22349SStefano 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);
2882df22349SStefano Zampini     ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr);
2892df22349SStefano Zampini   }
2902df22349SStefano Zampini   /* merge local matrices */
2912df22349SStefano Zampini   ii   = iptr;
2922df22349SStefano Zampini   jj   = jptr;
2932df22349SStefano Zampini   aa   = data;
2942df22349SStefano Zampini   *ii  = *(hdi++) + *(hoi++);
2952df22349SStefano Zampini   for (jd=0,jo=0,cum=0; *ii<nnz; cum++) {
2962df22349SStefano Zampini     PetscScalar *aold = aa;
2972df22349SStefano Zampini     PetscInt    *jold = jj,nc = jd+jo;
2982df22349SStefano Zampini     for (; jd<*hdi; jd++) { *jj++ = *hdj++;      *aa++ = *hdd++; }
2992df22349SStefano Zampini     for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
3002df22349SStefano Zampini     *(++ii) = *(hdi++) + *(hoi++);
3012df22349SStefano Zampini     ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr);
3022df22349SStefano Zampini   }
3032df22349SStefano Zampini   for (; cum<dr; cum++) *(++ii) = nnz;
3042df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
305a033916dSStefano Zampini     Mat_SeqAIJ* a;
306a033916dSStefano Zampini 
3072df22349SStefano Zampini     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr);
3082df22349SStefano Zampini     ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr);
309a033916dSStefano Zampini     /* hack SeqAIJ */
310a033916dSStefano Zampini     a          = (Mat_SeqAIJ*)(lA->data);
311a033916dSStefano Zampini     a->free_a  = PETSC_TRUE;
312a033916dSStefano Zampini     a->free_ij = PETSC_TRUE;
3132df22349SStefano Zampini     ierr = MatDestroy(&lA);CHKERRQ(ierr);
3142df22349SStefano Zampini   }
3152df22349SStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3162df22349SStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3172df22349SStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
3182df22349SStefano Zampini     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
3192df22349SStefano Zampini   }
3202df22349SStefano Zampini   PetscFunctionReturn(0);
3212df22349SStefano Zampini }
3222df22349SStefano Zampini 
32363c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
32463c07aadSStefano Zampini {
32563c07aadSStefano Zampini   Mat_HYPRE      *hB;
32663c07aadSStefano Zampini   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
32763c07aadSStefano Zampini   PetscErrorCode ierr;
32863c07aadSStefano Zampini 
32963c07aadSStefano Zampini   PetscFunctionBegin;
33063c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
33163c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
33263c07aadSStefano Zampini     /* always destroy the old matrix and create a new memory;
33363c07aadSStefano Zampini        hope this does not churn the memory too much. The problem
33463c07aadSStefano Zampini        is I do not know if it is possible to put the matrix back to
33563c07aadSStefano Zampini        its initial state so that we can directly copy the values
33663c07aadSStefano Zampini        the second time through. */
33763c07aadSStefano Zampini     hB = (Mat_HYPRE*)((*B)->data);
33863c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
33963c07aadSStefano Zampini   } else {
34063c07aadSStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
34163c07aadSStefano Zampini     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
34263c07aadSStefano Zampini     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
34363c07aadSStefano Zampini     hB   = (Mat_HYPRE*)((*B)->data);
34463c07aadSStefano Zampini   }
34563c07aadSStefano Zampini   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
34663c07aadSStefano Zampini   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
3474ec6421dSstefano_zampini   (*B)->preallocated = PETSC_TRUE;
34863c07aadSStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
34963c07aadSStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
35063c07aadSStefano Zampini   PetscFunctionReturn(0);
35163c07aadSStefano Zampini }
35263c07aadSStefano Zampini 
353ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
35463c07aadSStefano Zampini {
35563c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
35663c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
35763c07aadSStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
35863c07aadSStefano Zampini   MPI_Comm           comm;
35963c07aadSStefano Zampini   PetscScalar        *da,*oa,*aptr;
36063c07aadSStefano Zampini   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
36163c07aadSStefano Zampini   PetscInt           i,dnnz,onnz,m,n;
36258968eb6SStefano Zampini   HYPRE_Int          type;
36363c07aadSStefano Zampini   PetscMPIInt        size;
36463c07aadSStefano Zampini   PetscErrorCode     ierr;
36563c07aadSStefano Zampini 
36663c07aadSStefano Zampini   PetscFunctionBegin;
36763c07aadSStefano Zampini   comm = PetscObjectComm((PetscObject)A);
36863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
36963c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
37063c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
37163c07aadSStefano Zampini     PetscBool ismpiaij,isseqaij;
37263c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
3734099cc6bSBarry Smith     ierr = PetscObjectBaseTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
37463c07aadSStefano Zampini     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
37563c07aadSStefano Zampini   }
37663c07aadSStefano Zampini   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
37763c07aadSStefano Zampini 
37863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
37963c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(parcsr);
38063c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(parcsr);
38163c07aadSStefano Zampini   m     = hypre_CSRMatrixNumRows(hdiag);
38263c07aadSStefano Zampini   n     = hypre_CSRMatrixNumCols(hdiag);
38363c07aadSStefano Zampini   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
38463c07aadSStefano Zampini   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
385225daaf8SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
38663c07aadSStefano Zampini     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
38763c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
38863c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
389225daaf8SStefano Zampini   } else if (reuse == MAT_REUSE_MATRIX) {
39063c07aadSStefano Zampini     PetscInt  nr;
39163c07aadSStefano Zampini     PetscBool done;
39263c07aadSStefano Zampini     if (size > 1) {
39363c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
39463c07aadSStefano Zampini 
39563c07aadSStefano Zampini       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
39663c07aadSStefano 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);
39763c07aadSStefano 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);
39863c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
39963c07aadSStefano Zampini     } else {
40063c07aadSStefano Zampini       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
40163c07aadSStefano Zampini       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
40263c07aadSStefano 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);
40363c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
40463c07aadSStefano Zampini     }
405225daaf8SStefano Zampini   } else { /* MAT_INPLACE_MATRIX */
4067d968826Sstefano_zampini     dii = (PetscInt*)hypre_CSRMatrixI(hdiag);
4077d968826Sstefano_zampini     djj = (PetscInt*)hypre_CSRMatrixJ(hdiag);
408225daaf8SStefano Zampini     da  = hypre_CSRMatrixData(hdiag);
40963c07aadSStefano Zampini   }
41063c07aadSStefano Zampini   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
41163c07aadSStefano Zampini   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
41263c07aadSStefano Zampini   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
41363c07aadSStefano Zampini   iptr = djj;
41463c07aadSStefano Zampini   aptr = da;
41563c07aadSStefano Zampini   for (i=0; i<m; i++) {
41663c07aadSStefano Zampini     PetscInt nc = dii[i+1]-dii[i];
41763c07aadSStefano Zampini     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
41863c07aadSStefano Zampini     iptr += nc;
41963c07aadSStefano Zampini     aptr += nc;
42063c07aadSStefano Zampini   }
42163c07aadSStefano Zampini   if (size > 1) {
4227d968826Sstefano_zampini     HYPRE_Int *offdj,*coffd;
42363c07aadSStefano Zampini 
424225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
42563c07aadSStefano Zampini       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
42663c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
42763c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
428225daaf8SStefano Zampini     } else if (reuse == MAT_REUSE_MATRIX) {
42963c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
43063c07aadSStefano Zampini       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
43163c07aadSStefano Zampini       PetscBool  done;
43263c07aadSStefano Zampini 
43363c07aadSStefano Zampini       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
43463c07aadSStefano 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);
43563c07aadSStefano 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);
43663c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
437225daaf8SStefano Zampini     } else { /* MAT_INPLACE_MATRIX */
4387d968826Sstefano_zampini       oii = (PetscInt*)hypre_CSRMatrixI(hoffd);
4397d968826Sstefano_zampini       ojj = (PetscInt*)hypre_CSRMatrixJ(hoffd);
440225daaf8SStefano Zampini       oa  = hypre_CSRMatrixData(hoffd);
44163c07aadSStefano Zampini     }
44263c07aadSStefano Zampini     ierr  = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
44363c07aadSStefano Zampini     offdj = hypre_CSRMatrixJ(hoffd);
44463c07aadSStefano Zampini     coffd = hypre_ParCSRMatrixColMapOffd(parcsr);
44563c07aadSStefano Zampini     for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]];
44663c07aadSStefano Zampini     ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr);
44763c07aadSStefano Zampini     iptr = ojj;
44863c07aadSStefano Zampini     aptr = oa;
44963c07aadSStefano Zampini     for (i=0; i<m; i++) {
45063c07aadSStefano Zampini        PetscInt nc = oii[i+1]-oii[i];
45163c07aadSStefano Zampini        ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
45263c07aadSStefano Zampini        iptr += nc;
45363c07aadSStefano Zampini        aptr += nc;
45463c07aadSStefano Zampini     }
455225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
45663c07aadSStefano Zampini       Mat_MPIAIJ *b;
45763c07aadSStefano Zampini       Mat_SeqAIJ *d,*o;
458225daaf8SStefano Zampini 
45941073d71Sstefano_zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,B);CHKERRQ(ierr);
46063c07aadSStefano Zampini       /* hack MPIAIJ */
46163c07aadSStefano Zampini       b          = (Mat_MPIAIJ*)((*B)->data);
46263c07aadSStefano Zampini       d          = (Mat_SeqAIJ*)b->A->data;
46363c07aadSStefano Zampini       o          = (Mat_SeqAIJ*)b->B->data;
46463c07aadSStefano Zampini       d->free_a  = PETSC_TRUE;
46563c07aadSStefano Zampini       d->free_ij = PETSC_TRUE;
46663c07aadSStefano Zampini       o->free_a  = PETSC_TRUE;
46763c07aadSStefano Zampini       o->free_ij = PETSC_TRUE;
468225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
469225daaf8SStefano Zampini       Mat T;
47041073d71Sstefano_zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,&T);CHKERRQ(ierr);
471225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
472225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
473225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
474225daaf8SStefano Zampini       hypre_CSRMatrixI(hoffd)    = NULL;
475225daaf8SStefano Zampini       hypre_CSRMatrixJ(hoffd)    = NULL;
476225daaf8SStefano Zampini       hypre_CSRMatrixData(hoffd) = NULL;
477225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
47863c07aadSStefano Zampini     }
479225daaf8SStefano Zampini   } else {
480225daaf8SStefano Zampini     oii  = NULL;
481225daaf8SStefano Zampini     ojj  = NULL;
482225daaf8SStefano Zampini     oa   = NULL;
483225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
48463c07aadSStefano Zampini       Mat_SeqAIJ* b;
48563c07aadSStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
48663c07aadSStefano Zampini       /* hack SeqAIJ */
48763c07aadSStefano Zampini       b          = (Mat_SeqAIJ*)((*B)->data);
48863c07aadSStefano Zampini       b->free_a  = PETSC_TRUE;
48963c07aadSStefano Zampini       b->free_ij = PETSC_TRUE;
490225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
491225daaf8SStefano Zampini       Mat T;
492225daaf8SStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,&T);CHKERRQ(ierr);
493225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
494225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
495225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
496225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
49763c07aadSStefano Zampini     }
498225daaf8SStefano Zampini   }
499225daaf8SStefano Zampini 
500225daaf8SStefano Zampini   /* we have to use hypre_Tfree to free the arrays */
50163c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
502225daaf8SStefano Zampini     void *ptrs[6] = {dii,djj,da,oii,ojj,oa};
503225daaf8SStefano Zampini     const char *names[6] = {"_hypre_csr_dii",
504225daaf8SStefano Zampini                             "_hypre_csr_djj",
505225daaf8SStefano Zampini                             "_hypre_csr_da",
506225daaf8SStefano Zampini                             "_hypre_csr_oii",
507225daaf8SStefano Zampini                             "_hypre_csr_ojj",
508225daaf8SStefano Zampini                             "_hypre_csr_oa"};
509225daaf8SStefano Zampini     for (i=0; i<6; i++) {
510225daaf8SStefano Zampini       PetscContainer c;
511225daaf8SStefano Zampini 
512225daaf8SStefano Zampini       ierr = PetscContainerCreate(comm,&c);CHKERRQ(ierr);
513225daaf8SStefano Zampini       ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr);
514225daaf8SStefano Zampini       ierr = PetscContainerSetUserDestroy(c,hypre_array_destroy);CHKERRQ(ierr);
515225daaf8SStefano Zampini       ierr = PetscObjectCompose((PetscObject)(*B),names[i],(PetscObject)c);CHKERRQ(ierr);
516225daaf8SStefano Zampini       ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
517225daaf8SStefano Zampini     }
51863c07aadSStefano Zampini   }
51963c07aadSStefano Zampini   PetscFunctionReturn(0);
52063c07aadSStefano Zampini }
52163c07aadSStefano Zampini 
522613e5ff0Sstefano_zampini static PetscErrorCode MatAIJGetParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
523c1a070e6SStefano Zampini {
524613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *tA;
525c1a070e6SStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
526c1a070e6SStefano Zampini   Mat_SeqAIJ         *diag,*offd;
527c1a070e6SStefano Zampini   PetscInt           *garray,noffd,dnnz,onnz,*row_starts,*col_starts;
528c1a070e6SStefano Zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
529613e5ff0Sstefano_zampini   PetscBool          ismpiaij,isseqaij;
530c1a070e6SStefano Zampini   PetscErrorCode     ierr;
531c1a070e6SStefano Zampini 
532c1a070e6SStefano Zampini   PetscFunctionBegin;
533c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
5344099cc6bSBarry Smith   ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
535c1a070e6SStefano Zampini   if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type);
536c1a070e6SStefano Zampini   if (ismpiaij) {
537c1a070e6SStefano Zampini     Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data);
538c1a070e6SStefano Zampini 
539c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)a->A->data;
540c1a070e6SStefano Zampini     offd   = (Mat_SeqAIJ*)a->B->data;
541c1a070e6SStefano Zampini     garray = a->garray;
542c1a070e6SStefano Zampini     noffd  = a->B->cmap->N;
543c1a070e6SStefano Zampini     dnnz   = diag->nz;
544c1a070e6SStefano Zampini     onnz   = offd->nz;
545c1a070e6SStefano Zampini   } else {
546c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)A->data;
547c1a070e6SStefano Zampini     offd   = NULL;
548c1a070e6SStefano Zampini     garray = NULL;
549c1a070e6SStefano Zampini     noffd  = 0;
550c1a070e6SStefano Zampini     dnnz   = diag->nz;
551c1a070e6SStefano Zampini     onnz   = 0;
552c1a070e6SStefano Zampini   }
553225daaf8SStefano Zampini 
554c1a070e6SStefano Zampini   /* create a temporary ParCSR */
555c1a070e6SStefano Zampini   if (HYPRE_AssumedPartitionCheck()) {
556c1a070e6SStefano Zampini     PetscMPIInt myid;
557c1a070e6SStefano Zampini 
558c1a070e6SStefano Zampini     ierr       = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
559c1a070e6SStefano Zampini     row_starts = A->rmap->range + myid;
560c1a070e6SStefano Zampini     col_starts = A->cmap->range + myid;
561c1a070e6SStefano Zampini   } else {
562c1a070e6SStefano Zampini     row_starts = A->rmap->range;
563c1a070e6SStefano Zampini     col_starts = A->cmap->range;
564c1a070e6SStefano Zampini   }
5657d968826Sstefano_zampini   tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,(HYPRE_Int*)row_starts,(HYPRE_Int*)col_starts,noffd,dnnz,onnz);
566c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(tA,0);
567c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(tA,0);
568c1a070e6SStefano Zampini 
569225daaf8SStefano Zampini   /* set diagonal part */
570c1a070e6SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(tA);
5717d968826Sstefano_zampini   hypre_CSRMatrixI(hdiag)           = (HYPRE_Int*)diag->i;
5727d968826Sstefano_zampini   hypre_CSRMatrixJ(hdiag)           = (HYPRE_Int*)diag->j;
573c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = diag->a;
574c1a070e6SStefano Zampini   hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz;
575c1a070e6SStefano Zampini   hypre_CSRMatrixSetRownnz(hdiag);
576c1a070e6SStefano Zampini   hypre_CSRMatrixSetDataOwner(hdiag,0);
577c1a070e6SStefano Zampini 
578225daaf8SStefano Zampini   /* set offdiagonal part */
579c1a070e6SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(tA);
580c1a070e6SStefano Zampini   if (offd) {
5817d968826Sstefano_zampini     hypre_CSRMatrixI(hoffd)           = (HYPRE_Int*)offd->i;
5827d968826Sstefano_zampini     hypre_CSRMatrixJ(hoffd)           = (HYPRE_Int*)offd->j;
583c1a070e6SStefano Zampini     hypre_CSRMatrixData(hoffd)        = offd->a;
584c1a070e6SStefano Zampini     hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz;
585c1a070e6SStefano Zampini     hypre_CSRMatrixSetRownnz(hoffd);
586c1a070e6SStefano Zampini     hypre_CSRMatrixSetDataOwner(hoffd,0);
587c1a070e6SStefano Zampini     hypre_ParCSRMatrixSetNumNonzeros(tA);
5887d968826Sstefano_zampini     hypre_ParCSRMatrixColMapOffd(tA) = (HYPRE_Int*)garray;
589c1a070e6SStefano Zampini   }
590613e5ff0Sstefano_zampini   *hA = tA;
591613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
592613e5ff0Sstefano_zampini }
593c1a070e6SStefano Zampini 
594613e5ff0Sstefano_zampini static PetscErrorCode MatAIJRestoreParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
595613e5ff0Sstefano_zampini {
596613e5ff0Sstefano_zampini   hypre_CSRMatrix    *hdiag,*hoffd;
597c1a070e6SStefano Zampini 
598613e5ff0Sstefano_zampini   PetscFunctionBegin;
599613e5ff0Sstefano_zampini   hdiag = hypre_ParCSRMatrixDiag(*hA);
600613e5ff0Sstefano_zampini   hoffd = hypre_ParCSRMatrixOffd(*hA);
601c1a070e6SStefano Zampini   /* set pointers to NULL before destroying tA */
602c1a070e6SStefano Zampini   hypre_CSRMatrixI(hdiag)           = NULL;
603c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hdiag)           = NULL;
604c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = NULL;
605c1a070e6SStefano Zampini   hypre_CSRMatrixI(hoffd)           = NULL;
606c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hoffd)           = NULL;
607c1a070e6SStefano Zampini   hypre_CSRMatrixData(hoffd)        = NULL;
608613e5ff0Sstefano_zampini   hypre_ParCSRMatrixColMapOffd(*hA) = NULL;
609613e5ff0Sstefano_zampini   hypre_ParCSRMatrixDestroy(*hA);
610613e5ff0Sstefano_zampini   *hA = NULL;
611613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
612613e5ff0Sstefano_zampini }
613613e5ff0Sstefano_zampini 
614613e5ff0Sstefano_zampini /* calls RAP from BoomerAMG:
6153dad0653Sstefano_zampini    the resulting ParCSR will not own the column and row starts
6163dad0653Sstefano_zampini    It looks like we don't need to have the diagonal entries
6173dad0653Sstefano_zampini    ordered first in the rows of the diagonal part
6183dad0653Sstefano_zampini    for boomerAMGBuildCoarseOperator to work */
619a055b5aaSBarry Smith static PetscErrorCode MatHYPRE_ParCSR_RAP(hypre_ParCSRMatrix *hR, hypre_ParCSRMatrix *hA,hypre_ParCSRMatrix *hP, hypre_ParCSRMatrix **hRAP)
620613e5ff0Sstefano_zampini {
621613e5ff0Sstefano_zampini   PetscErrorCode ierr;
622613e5ff0Sstefano_zampini   HYPRE_Int      P_owns_col_starts,R_owns_row_starts;
623613e5ff0Sstefano_zampini 
624613e5ff0Sstefano_zampini   PetscFunctionBegin;
625613e5ff0Sstefano_zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(hP);
626613e5ff0Sstefano_zampini   R_owns_row_starts = hypre_ParCSRMatrixOwnsRowStarts(hR);
627613e5ff0Sstefano_zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(hR,hA,hP,hRAP));
628613e5ff0Sstefano_zampini   PetscStackCallStandard(hypre_ParCSRMatrixSetNumNonzeros,(*hRAP));
629613e5ff0Sstefano_zampini   /* hypre_BoomerAMGBuildCoarseOperator steals the col_starts from P and the row_starts from R */
630613e5ff0Sstefano_zampini   hypre_ParCSRMatrixSetRowStartsOwner(*hRAP,0);
631613e5ff0Sstefano_zampini   hypre_ParCSRMatrixSetColStartsOwner(*hRAP,0);
632613e5ff0Sstefano_zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(hP,1);
633613e5ff0Sstefano_zampini   if (R_owns_row_starts) hypre_ParCSRMatrixSetRowStartsOwner(hR,1);
634613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
635613e5ff0Sstefano_zampini }
636613e5ff0Sstefano_zampini 
6376f231fbdSstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat P,Mat C)
638613e5ff0Sstefano_zampini {
6396f231fbdSstefano_zampini   Mat                B;
640613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *hA,*hP,*hPtAP;
641613e5ff0Sstefano_zampini   PetscErrorCode     ierr;
642613e5ff0Sstefano_zampini 
643613e5ff0Sstefano_zampini   PetscFunctionBegin;
644613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
645613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(P,&hP);CHKERRQ(ierr);
646613e5ff0Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(hP,hA,hP,&hPtAP);CHKERRQ(ierr);
6476f231fbdSstefano_zampini   ierr = MatCreateFromParCSR(hPtAP,MATAIJ,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
6486f231fbdSstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
649613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
650613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(P,&hP);CHKERRQ(ierr);
6516f231fbdSstefano_zampini   PetscFunctionReturn(0);
6526f231fbdSstefano_zampini }
6536f231fbdSstefano_zampini 
6546f231fbdSstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat P,PetscReal fill,Mat *C)
6556f231fbdSstefano_zampini {
6566f231fbdSstefano_zampini   PetscErrorCode ierr;
657ab4d48faSStefano Zampini 
6586f231fbdSstefano_zampini   PetscFunctionBegin;
6596f231fbdSstefano_zampini   ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
6606f231fbdSstefano_zampini   ierr                   = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
6616f231fbdSstefano_zampini   (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_AIJ_wHYPRE;
662613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
663613e5ff0Sstefano_zampini }
664613e5ff0Sstefano_zampini 
6654cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_HYPRE(Mat A,Mat P,Mat C)
666613e5ff0Sstefano_zampini {
6674cc28894Sstefano_zampini   Mat                B;
6684cc28894Sstefano_zampini   Mat_HYPRE          *hP;
669613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *hA = NULL,*Pparcsr,*ptapparcsr;
670613e5ff0Sstefano_zampini   HYPRE_Int          type;
671613e5ff0Sstefano_zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
6724cc28894Sstefano_zampini   PetscBool          ishypre;
673613e5ff0Sstefano_zampini   PetscErrorCode     ierr;
674613e5ff0Sstefano_zampini 
675613e5ff0Sstefano_zampini   PetscFunctionBegin;
6764cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
6774cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(comm,PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
6784cc28894Sstefano_zampini   hP = (Mat_HYPRE*)P->data;
679613e5ff0Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
680613e5ff0Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
681613e5ff0Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
682613e5ff0Sstefano_zampini 
683613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
684613e5ff0Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,hA,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
685613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
686225daaf8SStefano Zampini 
6874cc28894Sstefano_zampini   /* create temporary matrix and merge to C */
6884cc28894Sstefano_zampini   ierr = MatCreateFromParCSR(ptapparcsr,((PetscObject)C)->type_name,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
6894cc28894Sstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
6904cc28894Sstefano_zampini   PetscFunctionReturn(0);
6914cc28894Sstefano_zampini }
6924cc28894Sstefano_zampini 
6934cc28894Sstefano_zampini static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
6944cc28894Sstefano_zampini {
6954cc28894Sstefano_zampini   PetscErrorCode ierr;
6964cc28894Sstefano_zampini 
6974cc28894Sstefano_zampini   PetscFunctionBegin;
6984cc28894Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
6994cc28894Sstefano_zampini     const char *deft = MATAIJ;
7004cc28894Sstefano_zampini     char       type[256];
7014cc28894Sstefano_zampini     PetscBool  flg;
7024cc28894Sstefano_zampini 
7034cc28894Sstefano_zampini     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
7044cc28894Sstefano_zampini     ierr = PetscOptionsFList("-matptap_hypre_outtype","Matrix type",NULL,MatList,deft,type,256,&flg);CHKERRQ(ierr);
7054cc28894Sstefano_zampini     ierr = PetscOptionsEnd();CHKERRQ(ierr);
7064cc28894Sstefano_zampini     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7074cc28894Sstefano_zampini     ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
7084cc28894Sstefano_zampini     if (flg) {
7094cc28894Sstefano_zampini       ierr = MatSetType(*C,type);CHKERRQ(ierr);
7104cc28894Sstefano_zampini     } else {
7114cc28894Sstefano_zampini       ierr = MatSetType(*C,deft);CHKERRQ(ierr);
7124cc28894Sstefano_zampini     }
7134cc28894Sstefano_zampini     (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_HYPRE;
7144cc28894Sstefano_zampini     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7154cc28894Sstefano_zampini   }
7164cc28894Sstefano_zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
7174cc28894Sstefano_zampini   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
718613e5ff0Sstefano_zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
719c1a070e6SStefano Zampini   PetscFunctionReturn(0);
720c1a070e6SStefano Zampini }
721c1a070e6SStefano Zampini 
7224cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_HYPRE_HYPRE(Mat A,Mat P,Mat C)
7234cc28894Sstefano_zampini {
7244cc28894Sstefano_zampini   Mat                B;
7254cc28894Sstefano_zampini   hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr;
7264cc28894Sstefano_zampini   Mat_HYPRE          *hA,*hP;
7274cc28894Sstefano_zampini   PetscBool          ishypre;
7284cc28894Sstefano_zampini   HYPRE_Int          type;
7294cc28894Sstefano_zampini   PetscErrorCode     ierr;
7304cc28894Sstefano_zampini 
7314cc28894Sstefano_zampini   PetscFunctionBegin;
7324cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
7334cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)P),PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
7344cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr);
7354cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE);
7364cc28894Sstefano_zampini   hA = (Mat_HYPRE*)A->data;
7374cc28894Sstefano_zampini   hP = (Mat_HYPRE*)P->data;
7384cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
7394cc28894Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
7404cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
7414cc28894Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
7424cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
7434cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
7444cc28894Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
7454cc28894Sstefano_zampini   ierr = MatCreateFromParCSR(ptapparcsr,MATHYPRE,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
7464cc28894Sstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
7474cc28894Sstefano_zampini   PetscFunctionReturn(0);
7484cc28894Sstefano_zampini }
7494cc28894Sstefano_zampini 
750cd8bc7baSStefano Zampini static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
751cd8bc7baSStefano Zampini {
752cd8bc7baSStefano Zampini   PetscErrorCode ierr;
753cd8bc7baSStefano Zampini 
754cd8bc7baSStefano Zampini   PetscFunctionBegin;
7554cc28894Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
7564cc28894Sstefano_zampini     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7574cc28894Sstefano_zampini     ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
7584cc28894Sstefano_zampini     ierr                   = MatSetType(*C,MATHYPRE);CHKERRQ(ierr);
7594cc28894Sstefano_zampini     (*C)->ops->ptapnumeric = MatPtAPNumeric_HYPRE_HYPRE;
7604cc28894Sstefano_zampini     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7614cc28894Sstefano_zampini   }
762613e5ff0Sstefano_zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
7634cc28894Sstefano_zampini   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
764613e5ff0Sstefano_zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
765cd8bc7baSStefano Zampini   PetscFunctionReturn(0);
766cd8bc7baSStefano Zampini }
767cd8bc7baSStefano Zampini 
768d501dc42Sstefano_zampini /* calls hypre_ParMatmul
769d501dc42Sstefano_zampini    hypre_ParMatMul uses hypre_ParMatrixCreate with the communicator of hA
7703dad0653Sstefano_zampini    hypre_ParMatrixCreate does not duplicate the communicator
7713dad0653Sstefano_zampini    It looks like we don't need to have the diagonal entries
7723dad0653Sstefano_zampini    ordered first in the rows of the diagonal part
7733dad0653Sstefano_zampini    for boomerAMGBuildCoarseOperator to work */
774d501dc42Sstefano_zampini static PetscErrorCode MatHYPRE_ParCSR_MatMatMult(hypre_ParCSRMatrix *hA, hypre_ParCSRMatrix *hB, hypre_ParCSRMatrix **hAB)
775d501dc42Sstefano_zampini {
776d501dc42Sstefano_zampini   PetscFunctionBegin;
777d501dc42Sstefano_zampini   PetscStackPush("hypre_ParMatmul");
778d501dc42Sstefano_zampini   *hAB = hypre_ParMatmul(hA,hB);
779d501dc42Sstefano_zampini   PetscStackPop;
780d501dc42Sstefano_zampini   PetscFunctionReturn(0);
781d501dc42Sstefano_zampini }
782d501dc42Sstefano_zampini 
7835e5acdf2Sstefano_zampini static PetscErrorCode MatMatMultNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C)
7845e5acdf2Sstefano_zampini {
7855e5acdf2Sstefano_zampini   Mat                D;
786d501dc42Sstefano_zampini   hypre_ParCSRMatrix *hA,*hB,*hAB = NULL;
7875e5acdf2Sstefano_zampini   PetscErrorCode     ierr;
7885e5acdf2Sstefano_zampini 
7895e5acdf2Sstefano_zampini   PetscFunctionBegin;
7905e5acdf2Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
7915e5acdf2Sstefano_zampini   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
792d501dc42Sstefano_zampini   ierr = MatHYPRE_ParCSR_MatMatMult(hA,hB,&hAB);CHKERRQ(ierr);
7935e5acdf2Sstefano_zampini   ierr = MatCreateFromParCSR(hAB,MATAIJ,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
7945e5acdf2Sstefano_zampini   ierr = MatHeaderMerge(C,&D);CHKERRQ(ierr);
7955e5acdf2Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
7965e5acdf2Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
7975e5acdf2Sstefano_zampini   PetscFunctionReturn(0);
7985e5acdf2Sstefano_zampini }
7995e5acdf2Sstefano_zampini 
8005e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat B,PetscReal fill,Mat *C)
8015e5acdf2Sstefano_zampini {
8025e5acdf2Sstefano_zampini   PetscErrorCode ierr;
80320e1dc0dSstefano_zampini 
8045e5acdf2Sstefano_zampini   PetscFunctionBegin;
8055e5acdf2Sstefano_zampini   ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
8065e5acdf2Sstefano_zampini   ierr                      = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
8075e5acdf2Sstefano_zampini   (*C)->ops->matmultnumeric = MatMatMultNumeric_AIJ_AIJ_wHYPRE;
8085e5acdf2Sstefano_zampini   PetscFunctionReturn(0);
8095e5acdf2Sstefano_zampini }
8105e5acdf2Sstefano_zampini 
811d501dc42Sstefano_zampini static PetscErrorCode MatMatMultNumeric_HYPRE_HYPRE(Mat A,Mat B,Mat C)
812d501dc42Sstefano_zampini {
813d501dc42Sstefano_zampini   Mat                D;
814d501dc42Sstefano_zampini   hypre_ParCSRMatrix *Aparcsr,*Bparcsr,*ABparcsr = NULL;
815d501dc42Sstefano_zampini   Mat_HYPRE          *hA,*hB;
816d501dc42Sstefano_zampini   PetscBool          ishypre;
817d501dc42Sstefano_zampini   HYPRE_Int          type;
818d501dc42Sstefano_zampini   PetscErrorCode     ierr;
819d501dc42Sstefano_zampini 
820d501dc42Sstefano_zampini   PetscFunctionBegin;
821d501dc42Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)B,MATHYPRE,&ishypre);CHKERRQ(ierr);
822d501dc42Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_USER,"B should be of type %s",MATHYPRE);
823d501dc42Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr);
824d501dc42Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE);
825d501dc42Sstefano_zampini   hA = (Mat_HYPRE*)A->data;
826d501dc42Sstefano_zampini   hB = (Mat_HYPRE*)B->data;
827d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
828d501dc42Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
829d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hB->ij,&type));
830d501dc42Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
831d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
832d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hB->ij,(void**)&Bparcsr));
833d501dc42Sstefano_zampini   ierr = MatHYPRE_ParCSR_MatMatMult(Aparcsr,Bparcsr,&ABparcsr);CHKERRQ(ierr);
834d501dc42Sstefano_zampini   ierr = MatCreateFromParCSR(ABparcsr,MATHYPRE,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
835d501dc42Sstefano_zampini   /* need to use HeaderReplace because HeaderMerge messes up with the communicator */
836d501dc42Sstefano_zampini   ierr = MatHeaderReplace(C,&D);CHKERRQ(ierr);
837d501dc42Sstefano_zampini   C->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE;
838d501dc42Sstefano_zampini   PetscFunctionReturn(0);
839d501dc42Sstefano_zampini }
840d501dc42Sstefano_zampini 
841d501dc42Sstefano_zampini static PetscErrorCode MatMatMult_HYPRE_HYPRE(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
842d501dc42Sstefano_zampini {
843d501dc42Sstefano_zampini   PetscErrorCode ierr;
844d501dc42Sstefano_zampini 
845d501dc42Sstefano_zampini   PetscFunctionBegin;
846d501dc42Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
847d501dc42Sstefano_zampini     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
848d501dc42Sstefano_zampini     ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
849d501dc42Sstefano_zampini     ierr                      = MatSetType(*C,MATHYPRE);CHKERRQ(ierr);
850d501dc42Sstefano_zampini     (*C)->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE;
851d501dc42Sstefano_zampini     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
852d501dc42Sstefano_zampini   }
853d501dc42Sstefano_zampini   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
854d501dc42Sstefano_zampini   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
855d501dc42Sstefano_zampini   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
856d501dc42Sstefano_zampini   PetscFunctionReturn(0);
857d501dc42Sstefano_zampini }
858d501dc42Sstefano_zampini 
8593dad0653Sstefano_zampini PETSC_INTERN PetscErrorCode MatTransposeMatMatMultNumeric_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,Mat D)
86020e1dc0dSstefano_zampini {
86120e1dc0dSstefano_zampini   Mat                E;
86220e1dc0dSstefano_zampini   hypre_ParCSRMatrix *hA,*hB,*hC,*hABC;
86320e1dc0dSstefano_zampini   PetscErrorCode     ierr;
86420e1dc0dSstefano_zampini 
86520e1dc0dSstefano_zampini   PetscFunctionBegin;
86620e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
86720e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
86820e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(C,&hC);CHKERRQ(ierr);
86920e1dc0dSstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(hA,hB,hC,&hABC);CHKERRQ(ierr);
87020e1dc0dSstefano_zampini   ierr = MatCreateFromParCSR(hABC,MATAIJ,PETSC_OWN_POINTER,&E);CHKERRQ(ierr);
87120e1dc0dSstefano_zampini   ierr = MatHeaderMerge(D,&E);CHKERRQ(ierr);
87220e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
87320e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
87420e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(C,&hC);CHKERRQ(ierr);
87520e1dc0dSstefano_zampini   PetscFunctionReturn(0);
87620e1dc0dSstefano_zampini }
87720e1dc0dSstefano_zampini 
8783dad0653Sstefano_zampini PETSC_INTERN PetscErrorCode MatTransposeMatMatMultSymbolic_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,PetscReal fill,Mat *D)
87920e1dc0dSstefano_zampini {
88020e1dc0dSstefano_zampini   PetscErrorCode ierr;
88120e1dc0dSstefano_zampini 
88220e1dc0dSstefano_zampini   PetscFunctionBegin;
88320e1dc0dSstefano_zampini   ierr = MatCreate(PetscObjectComm((PetscObject)A),D);CHKERRQ(ierr);
88420e1dc0dSstefano_zampini   ierr = MatSetType(*D,MATAIJ);CHKERRQ(ierr);
88520e1dc0dSstefano_zampini   PetscFunctionReturn(0);
88620e1dc0dSstefano_zampini }
88720e1dc0dSstefano_zampini 
888ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
88963c07aadSStefano Zampini {
89063c07aadSStefano Zampini   PetscErrorCode ierr;
89163c07aadSStefano Zampini 
89263c07aadSStefano Zampini   PetscFunctionBegin;
89363c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
89463c07aadSStefano Zampini   PetscFunctionReturn(0);
89563c07aadSStefano Zampini }
89663c07aadSStefano Zampini 
897ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
89863c07aadSStefano Zampini {
89963c07aadSStefano Zampini   PetscErrorCode ierr;
90063c07aadSStefano Zampini 
90163c07aadSStefano Zampini   PetscFunctionBegin;
90263c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
90363c07aadSStefano Zampini   PetscFunctionReturn(0);
90463c07aadSStefano Zampini }
90563c07aadSStefano Zampini 
906ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
90763c07aadSStefano Zampini {
90863c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
90963c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
91063c07aadSStefano Zampini   hypre_ParVector    *hx,*hy;
91163c07aadSStefano Zampini   PetscScalar        *ax,*ay,*sax,*say;
91263c07aadSStefano Zampini   PetscErrorCode     ierr;
91363c07aadSStefano Zampini 
91463c07aadSStefano Zampini   PetscFunctionBegin;
91563c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
91663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
91763c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
91863c07aadSStefano Zampini   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
91963c07aadSStefano Zampini   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
92063c07aadSStefano Zampini   if (trans) {
92158968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ay,say);
92258968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax);
92363c07aadSStefano Zampini     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
92458968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,say,ay);
92558968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax);
92663c07aadSStefano Zampini   } else {
92758968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax);
92858968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ay,say);
92963c07aadSStefano Zampini     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
93058968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax);
93158968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,say,ay);
93263c07aadSStefano Zampini   }
93363c07aadSStefano Zampini   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
93463c07aadSStefano Zampini   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
93563c07aadSStefano Zampini   PetscFunctionReturn(0);
93663c07aadSStefano Zampini }
93763c07aadSStefano Zampini 
938ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A)
93963c07aadSStefano Zampini {
94063c07aadSStefano Zampini   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
94163c07aadSStefano Zampini   PetscErrorCode ierr;
94263c07aadSStefano Zampini 
94363c07aadSStefano Zampini   PetscFunctionBegin;
94463c07aadSStefano Zampini   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
94563c07aadSStefano Zampini   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
946978814f1SStefano Zampini   if (hA->ij) {
947978814f1SStefano Zampini     if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL;
948978814f1SStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
949978814f1SStefano Zampini   }
95063c07aadSStefano Zampini   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); }
951c69f721fSFande Kong 
952c69f721fSFande Kong   ierr = MatStashDestroy_Private(&A->stash);CHKERRQ(ierr);
953c69f721fSFande Kong 
954c69f721fSFande Kong   ierr = PetscFree(hA->array);CHKERRQ(ierr);
955c69f721fSFande Kong 
95663c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
9572df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr);
958c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr);
959c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr);
960d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPRESetPreallocation_C",NULL);CHKERRQ(ierr);
961dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPREGetParCSR_C",NULL);CHKERRQ(ierr);
96263c07aadSStefano Zampini   ierr = PetscFree(A->data);CHKERRQ(ierr);
96363c07aadSStefano Zampini   PetscFunctionReturn(0);
96463c07aadSStefano Zampini }
96563c07aadSStefano Zampini 
966ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A)
96763c07aadSStefano Zampini {
9684ec6421dSstefano_zampini   PetscErrorCode ierr;
9694ec6421dSstefano_zampini 
9704ec6421dSstefano_zampini   PetscFunctionBegin;
9714ec6421dSstefano_zampini   ierr = MatHYPRESetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr);
9724ec6421dSstefano_zampini   PetscFunctionReturn(0);
9734ec6421dSstefano_zampini }
9744ec6421dSstefano_zampini 
975c69f721fSFande Kong 
9764ec6421dSstefano_zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
9774ec6421dSstefano_zampini {
97863c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
97963c07aadSStefano Zampini   Vec                x,b;
980c69f721fSFande Kong   PetscMPIInt        n;
981c69f721fSFande Kong   PetscInt           i,j,rstart,ncols,flg;
982c69f721fSFande Kong   PetscInt           *row,*col;
983c69f721fSFande Kong   PetscScalar        *val;
98463c07aadSStefano Zampini   PetscErrorCode     ierr;
98563c07aadSStefano Zampini 
98663c07aadSStefano Zampini   PetscFunctionBegin;
9874ec6421dSstefano_zampini   if (mode == MAT_FLUSH_ASSEMBLY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_FLUSH_ASSEMBLY currently not supported with MATHYPRE");
988c69f721fSFande Kong 
989c69f721fSFande Kong   if (!A->nooffprocentries) {
990c69f721fSFande Kong     while (1) {
991c69f721fSFande Kong       ierr = MatStashScatterGetMesg_Private(&A->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
992c69f721fSFande Kong       if (!flg) break;
993c69f721fSFande Kong 
994c69f721fSFande Kong       for (i=0; i<n; ) {
995c69f721fSFande Kong         /* Now identify the consecutive vals belonging to the same row */
996c69f721fSFande Kong         for (j=i,rstart=row[j]; j<n; j++) {
997c69f721fSFande Kong           if (row[j] != rstart) break;
998c69f721fSFande Kong         }
999c69f721fSFande Kong         if (j < n) ncols = j-i;
1000c69f721fSFande Kong         else       ncols = n-i;
1001c69f721fSFande Kong         /* Now assemble all these values with a single function call */
1002c69f721fSFande Kong         ierr = MatSetValues_HYPRE(A,1,row+i,ncols,col+i,val+i,A->insertmode);CHKERRQ(ierr);
1003c69f721fSFande Kong 
1004c69f721fSFande Kong         i = j;
1005c69f721fSFande Kong       }
1006c69f721fSFande Kong     }
1007c69f721fSFande Kong     ierr = MatStashScatterEnd_Private(&A->stash);CHKERRQ(ierr);
1008c69f721fSFande Kong   }
1009c69f721fSFande Kong 
10104ec6421dSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
101163c07aadSStefano Zampini   if (hA->x) PetscFunctionReturn(0);
101263c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
101363c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
101463c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
101563c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
101663c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
101763c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
101863c07aadSStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
101963c07aadSStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
102063c07aadSStefano Zampini   PetscFunctionReturn(0);
102163c07aadSStefano Zampini }
102263c07aadSStefano Zampini 
1023c69f721fSFande Kong static PetscErrorCode MatGetArray_HYPRE(Mat A, PetscInt size, void **array)
1024c69f721fSFande Kong {
1025c69f721fSFande Kong   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
1026c69f721fSFande Kong   PetscErrorCode     ierr;
1027c69f721fSFande Kong 
1028c69f721fSFande Kong   PetscFunctionBegin;
1029c69f721fSFande Kong   if (!hA->available) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Temporary space is in use");
1030c69f721fSFande Kong 
1031c69f721fSFande Kong   if (hA->size >= size) *array = hA->array;
1032c69f721fSFande Kong   else {
1033c69f721fSFande Kong     ierr = PetscFree(hA->array);CHKERRQ(ierr);
1034c69f721fSFande Kong     hA->size = size;
1035c69f721fSFande Kong     ierr = PetscMalloc(hA->size,&hA->array);CHKERRQ(ierr);
1036c69f721fSFande Kong     *array = hA->array;
1037c69f721fSFande Kong   }
1038c69f721fSFande Kong 
1039c69f721fSFande Kong   hA->available = PETSC_FALSE;
1040c69f721fSFande Kong   PetscFunctionReturn(0);
1041c69f721fSFande Kong }
1042c69f721fSFande Kong 
1043c69f721fSFande Kong static PetscErrorCode MatRestoreArray_HYPRE(Mat A, PetscInt size, void **array)
1044c69f721fSFande Kong {
1045c69f721fSFande Kong   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
1046c69f721fSFande Kong 
1047c69f721fSFande Kong   PetscFunctionBegin;
1048c69f721fSFande Kong   *array = NULL;
1049c69f721fSFande Kong   hA->available = PETSC_TRUE;
1050c69f721fSFande Kong   PetscFunctionReturn(0);
1051c69f721fSFande Kong }
1052c69f721fSFande Kong 
1053d975228cSstefano_zampini 
1054d975228cSstefano_zampini PetscErrorCode MatSetValues_HYPRE(Mat A, PetscInt nr, const PetscInt rows[], PetscInt nc, const PetscInt cols[], const PetscScalar v[], InsertMode ins)
1055d975228cSstefano_zampini {
1056d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
1057d975228cSstefano_zampini   PetscScalar        *vals = (PetscScalar *)v;
1058c69f721fSFande Kong   PetscScalar        *sscr;
1059c69f721fSFande Kong   PetscInt           *cscr[2];
1060c69f721fSFande Kong   PetscInt           i,nzc;
1061*08defe43SFande Kong   void               *array = NULL;
1062d975228cSstefano_zampini   PetscErrorCode     ierr;
1063d975228cSstefano_zampini 
1064d975228cSstefano_zampini   PetscFunctionBegin;
1065c69f721fSFande Kong   ierr = MatGetArray_HYPRE(A,sizeof(PetscInt)*(2*nc)+sizeof(PetscScalar)*nc*nr,&array);CHKERRQ(ierr);
1066*08defe43SFande Kong   ierr = PetscMemzero(array,sizeof(PetscInt)*(2*nc)+sizeof(PetscScalar)*nc*nr);CHKERRQ(ierr);
1067c69f721fSFande Kong   cscr[0] = (PetscInt*)array;
1068c69f721fSFande Kong   cscr[1] = ((PetscInt*)array)+nc;
1069c69f721fSFande Kong   sscr = (PetscScalar*)(((PetscInt*)array)+nc*2);
1070d975228cSstefano_zampini   for (i=0,nzc=0;i<nc;i++) {
1071d975228cSstefano_zampini     if (cols[i] >= 0) {
1072d975228cSstefano_zampini       cscr[0][nzc  ] = cols[i];
1073d975228cSstefano_zampini       cscr[1][nzc++] = i;
1074d975228cSstefano_zampini     }
1075d975228cSstefano_zampini   }
1076c69f721fSFande Kong   if (!nzc) {
1077c69f721fSFande Kong     ierr = MatRestoreArray_HYPRE(A,sizeof(PetscInt)*(2*nc)+sizeof(PetscScalar)*nc*nr,&array);CHKERRQ(ierr);
1078c69f721fSFande Kong     PetscFunctionReturn(0);
1079c69f721fSFande Kong   }
1080d975228cSstefano_zampini 
1081d975228cSstefano_zampini   if (ins == ADD_VALUES) {
1082d975228cSstefano_zampini     for (i=0;i<nr;i++) {
1083e3977e59Sstefano_zampini       if (rows[i] >= 0 && nzc) {
1084d975228cSstefano_zampini         PetscInt j;
1085d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
1086*08defe43SFande Kong         PetscStackCallStandard(HYPRE_IJMatrixAddToValues,(hA->ij,1,(HYPRE_Int*)&nzc,(HYPRE_Int*)(rows+i),(HYPRE_Int*)cscr[0],sscr));
1087d975228cSstefano_zampini       }
1088d975228cSstefano_zampini       vals += nc;
1089d975228cSstefano_zampini     }
1090d975228cSstefano_zampini   } else { /* INSERT_VALUES */
1091c69f721fSFande Kong 
1092d975228cSstefano_zampini     PetscInt rst,ren;
1093d975228cSstefano_zampini     ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
1094c69f721fSFande Kong 
1095d975228cSstefano_zampini     for (i=0;i<nr;i++) {
1096e3977e59Sstefano_zampini       if (rows[i] >= 0 && nzc) {
1097d975228cSstefano_zampini         PetscInt j;
1098d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
1099c69f721fSFande Kong         /* nonlocal values */
1100c69f721fSFande Kong         if (rows[i] < rst || rows[i] >= ren) { ierr = MatStashValuesRow_Private(&A->stash,rows[i],nzc,cscr[0],sscr, PETSC_FALSE);CHKERRQ(ierr); }
1101c69f721fSFande Kong         /* local values */
1102*08defe43SFande Kong         else PetscStackCallStandard(HYPRE_IJMatrixSetValues,(hA->ij,1,(HYPRE_Int*)&nzc,(HYPRE_Int*)(rows+i),(HYPRE_Int*)cscr[0],sscr));
1103d975228cSstefano_zampini       }
1104d975228cSstefano_zampini       vals += nc;
1105d975228cSstefano_zampini     }
1106d975228cSstefano_zampini   }
1107c69f721fSFande Kong 
1108c69f721fSFande Kong   ierr = MatRestoreArray_HYPRE(A,sizeof(PetscInt)*(2*nc)+sizeof(PetscScalar)*nc*nr,&array);CHKERRQ(ierr);
1109d975228cSstefano_zampini   PetscFunctionReturn(0);
1110d975228cSstefano_zampini }
1111d975228cSstefano_zampini 
1112d975228cSstefano_zampini static PetscErrorCode MatHYPRESetPreallocation_HYPRE(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1113d975228cSstefano_zampini {
1114d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
11157d968826Sstefano_zampini   HYPRE_Int          *hdnnz,*honnz;
111606a29025Sstefano_zampini   PetscInt           i,rs,re,cs,ce,bs;
1117d975228cSstefano_zampini   PetscMPIInt        size;
1118d975228cSstefano_zampini   PetscErrorCode     ierr;
1119d975228cSstefano_zampini 
1120d975228cSstefano_zampini   PetscFunctionBegin;
112106a29025Sstefano_zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1122d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1123d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1124d975228cSstefano_zampini   rs   = A->rmap->rstart;
1125d975228cSstefano_zampini   re   = A->rmap->rend;
1126d975228cSstefano_zampini   cs   = A->cmap->rstart;
1127d975228cSstefano_zampini   ce   = A->cmap->rend;
1128d975228cSstefano_zampini   if (!hA->ij) {
1129d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rs,re-1,cs,ce-1,&hA->ij));
1130d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1131d975228cSstefano_zampini   } else {
11327d968826Sstefano_zampini     HYPRE_Int hrs,hre,hcs,hce;
1133d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixGetLocalRange,(hA->ij,&hrs,&hre,&hcs,&hce));
1134d975228cSstefano_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);
1135d975228cSstefano_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);
1136d975228cSstefano_zampini   }
1137d975228cSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1138d975228cSstefano_zampini 
113906a29025Sstefano_zampini   if (dnz == PETSC_DEFAULT || dnz == PETSC_DECIDE) dnz = 10*bs;
114006a29025Sstefano_zampini   if (onz == PETSC_DEFAULT || onz == PETSC_DECIDE) onz = 10*bs;
114106a29025Sstefano_zampini 
1142d975228cSstefano_zampini   if (!dnnz) {
1143d975228cSstefano_zampini     ierr = PetscMalloc1(A->rmap->n,&hdnnz);CHKERRQ(ierr);
1144d975228cSstefano_zampini     for (i=0;i<A->rmap->n;i++) hdnnz[i] = dnz;
1145d975228cSstefano_zampini   } else {
11467d968826Sstefano_zampini     hdnnz = (HYPRE_Int*)dnnz;
1147d975228cSstefano_zampini   }
1148d975228cSstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1149d975228cSstefano_zampini   if (size > 1) {
1150d975228cSstefano_zampini     if (!onnz) {
1151d975228cSstefano_zampini       ierr = PetscMalloc1(A->rmap->n,&honnz);CHKERRQ(ierr);
1152d975228cSstefano_zampini       for (i=0;i<A->rmap->n;i++) honnz[i] = onz;
1153d975228cSstefano_zampini     } else {
11547d968826Sstefano_zampini       honnz = (HYPRE_Int*)onnz;
1155d975228cSstefano_zampini     }
1156d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(hA->ij,hdnnz,honnz));
1157d975228cSstefano_zampini   } else {
1158d975228cSstefano_zampini     honnz = NULL;
1159d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetRowSizes,(hA->ij,hdnnz));
1160d975228cSstefano_zampini   }
1161d975228cSstefano_zampini   if (!dnnz) {
1162d975228cSstefano_zampini     ierr = PetscFree(hdnnz);CHKERRQ(ierr);
1163d975228cSstefano_zampini   }
1164d975228cSstefano_zampini   if (!onnz && honnz) {
1165d975228cSstefano_zampini     ierr = PetscFree(honnz);CHKERRQ(ierr);
1166d975228cSstefano_zampini   }
116706a29025Sstefano_zampini   A->preallocated = PETSC_TRUE;
1168d975228cSstefano_zampini 
1169d975228cSstefano_zampini   /* SetDiagOffdSizes sets hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0 */
1170d975228cSstefano_zampini   {
1171d975228cSstefano_zampini     hypre_AuxParCSRMatrix *aux_matrix;
1172d975228cSstefano_zampini     aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij);
1173d975228cSstefano_zampini     hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1;
1174d975228cSstefano_zampini   }
1175d975228cSstefano_zampini   PetscFunctionReturn(0);
1176d975228cSstefano_zampini }
1177d975228cSstefano_zampini 
1178d975228cSstefano_zampini /*@C
1179d975228cSstefano_zampini    MatHYPRESetPreallocation - Preallocates memory for a sparse parallel matrix in HYPRE IJ format
1180d975228cSstefano_zampini 
1181d975228cSstefano_zampini    Collective on Mat
1182d975228cSstefano_zampini 
1183d975228cSstefano_zampini    Input Parameters:
1184d975228cSstefano_zampini +  A - the matrix
1185d975228cSstefano_zampini .  dnz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1186d975228cSstefano_zampini           (same value is used for all local rows)
1187d975228cSstefano_zampini .  dnnz - array containing the number of nonzeros in the various rows of the
1188d975228cSstefano_zampini           DIAGONAL portion of the local submatrix (possibly different for each row)
1189d975228cSstefano_zampini           or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure.
1190d975228cSstefano_zampini           The size of this array is equal to the number of local rows, i.e 'm'.
1191d975228cSstefano_zampini           For matrices that will be factored, you must leave room for (and set)
1192d975228cSstefano_zampini           the diagonal entry even if it is zero.
1193d975228cSstefano_zampini .  onz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1194d975228cSstefano_zampini           submatrix (same value is used for all local rows).
1195d975228cSstefano_zampini -  onnz - array containing the number of nonzeros in the various rows of the
1196d975228cSstefano_zampini           OFF-DIAGONAL portion of the local submatrix (possibly different for
1197d975228cSstefano_zampini           each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero
1198d975228cSstefano_zampini           structure. The size of this array is equal to the number
1199d975228cSstefano_zampini           of local rows, i.e 'm'.
1200d975228cSstefano_zampini 
1201d975228cSstefano_zampini    Notes: If the *nnz parameter is given then the *nz parameter is ignored; for sequential matrices, onz and onnz are ignored.
1202d975228cSstefano_zampini 
1203d975228cSstefano_zampini    Level: intermediate
1204d975228cSstefano_zampini 
1205d975228cSstefano_zampini .keywords: matrix, aij, compressed row, sparse, parallel
1206d975228cSstefano_zampini 
1207d975228cSstefano_zampini .seealso: MatCreate(), MatMPIAIJSetPreallocation, MATHYPRE
1208d975228cSstefano_zampini @*/
1209d975228cSstefano_zampini PetscErrorCode MatHYPRESetPreallocation(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1210d975228cSstefano_zampini {
1211d975228cSstefano_zampini   PetscErrorCode ierr;
1212d975228cSstefano_zampini 
1213d975228cSstefano_zampini   PetscFunctionBegin;
1214d975228cSstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1215d975228cSstefano_zampini   PetscValidType(A,1);
1216d975228cSstefano_zampini   ierr = PetscTryMethod(A,"MatHYPRESetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(A,dnz,dnnz,onz,onnz));CHKERRQ(ierr);
1217d975228cSstefano_zampini   PetscFunctionReturn(0);
1218d975228cSstefano_zampini }
1219d975228cSstefano_zampini 
1220225daaf8SStefano Zampini /*
1221225daaf8SStefano Zampini    MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix
1222225daaf8SStefano Zampini 
1223225daaf8SStefano Zampini    Collective
1224225daaf8SStefano Zampini 
1225225daaf8SStefano Zampini    Input Parameters:
1226225daaf8SStefano Zampini +  vparcsr  - the pointer to the hypre_ParCSRMatrix
1227bb4689ddSStefano Zampini .  mtype    - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported.
1228225daaf8SStefano Zampini -  copymode - PETSc copying options
1229225daaf8SStefano Zampini 
1230225daaf8SStefano Zampini    Output Parameter:
1231225daaf8SStefano Zampini .  A  - the matrix
1232225daaf8SStefano Zampini 
1233225daaf8SStefano Zampini    Level: intermediate
1234225daaf8SStefano Zampini 
1235225daaf8SStefano Zampini .seealso: MatHYPRE, PetscCopyMode
1236225daaf8SStefano Zampini */
1237dd9c0a25Sstefano_zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(hypre_ParCSRMatrix *vparcsr, MatType mtype, PetscCopyMode copymode, Mat* A)
1238978814f1SStefano Zampini {
1239225daaf8SStefano Zampini   Mat                   T;
1240978814f1SStefano Zampini   Mat_HYPRE             *hA;
124163c07aadSStefano Zampini   hypre_ParCSRMatrix    *parcsr;
1242978814f1SStefano Zampini   MPI_Comm              comm;
1243978814f1SStefano Zampini   PetscInt              rstart,rend,cstart,cend,M,N;
1244bb4689ddSStefano Zampini   PetscBool             isseqaij,ismpiaij,isaij,ishyp,isis;
1245978814f1SStefano Zampini   PetscErrorCode        ierr;
1246978814f1SStefano Zampini 
1247978814f1SStefano Zampini   PetscFunctionBegin;
1248978814f1SStefano Zampini   parcsr = (hypre_ParCSRMatrix *)vparcsr;
1249978814f1SStefano Zampini   comm   = hypre_ParCSRMatrixComm(parcsr);
1250225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
1251225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
1252225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr);
1253225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr);
12548cfe8d00SStefano Zampini   ierr   = PetscStrcmp(mtype,MATIS,&isis);CHKERRQ(ierr);
1255225daaf8SStefano Zampini   isaij  = (PetscBool)(isseqaij || ismpiaij || isaij);
1256bb4689ddSStefano 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);
1257225daaf8SStefano Zampini   if (ishyp && copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES");
1258978814f1SStefano Zampini 
1259978814f1SStefano Zampini   /* access ParCSRMatrix */
1260978814f1SStefano Zampini   rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr);
1261978814f1SStefano Zampini   rend   = hypre_ParCSRMatrixLastRowIndex(parcsr);
1262978814f1SStefano Zampini   cstart = hypre_ParCSRMatrixFirstColDiag(parcsr);
1263978814f1SStefano Zampini   cend   = hypre_ParCSRMatrixLastColDiag(parcsr);
1264978814f1SStefano Zampini   M      = hypre_ParCSRMatrixGlobalNumRows(parcsr);
1265978814f1SStefano Zampini   N      = hypre_ParCSRMatrixGlobalNumCols(parcsr);
1266978814f1SStefano Zampini 
1267fa92c42cSstefano_zampini   /* fix for empty local rows/columns */
1268fa92c42cSstefano_zampini   if (rend < rstart) rend = rstart;
1269fa92c42cSstefano_zampini   if (cend < cstart) cend = cstart;
1270fa92c42cSstefano_zampini 
1271978814f1SStefano Zampini   /* create PETSc matrix with MatHYPRE */
1272225daaf8SStefano Zampini   ierr = MatCreate(comm,&T);CHKERRQ(ierr);
1273225daaf8SStefano Zampini   ierr = MatSetSizes(T,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr);
1274225daaf8SStefano Zampini   ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr);
1275225daaf8SStefano Zampini   hA   = (Mat_HYPRE*)(T->data);
1276978814f1SStefano Zampini 
1277978814f1SStefano Zampini   /* create HYPRE_IJMatrix */
1278978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
1279978814f1SStefano Zampini 
1280978814f1SStefano Zampini   /* set ParCSR object */
1281978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1282978814f1SStefano Zampini   hypre_IJMatrixObject(hA->ij) = parcsr;
12834ec6421dSstefano_zampini   T->preallocated = PETSC_TRUE;
1284978814f1SStefano Zampini 
1285978814f1SStefano Zampini   /* set assembled flag */
1286978814f1SStefano Zampini   hypre_IJMatrixAssembleFlag(hA->ij) = 1;
1287978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1288225daaf8SStefano Zampini   if (ishyp) {
12896d2a658fSstefano_zampini     PetscMPIInt myid = 0;
12906d2a658fSstefano_zampini 
12916d2a658fSstefano_zampini     /* make sure we always have row_starts and col_starts available */
12926d2a658fSstefano_zampini     if (HYPRE_AssumedPartitionCheck()) {
12936d2a658fSstefano_zampini       ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
12946d2a658fSstefano_zampini     }
12956d2a658fSstefano_zampini     if (!hypre_ParCSRMatrixOwnsColStarts(parcsr)) {
12966d2a658fSstefano_zampini       PetscLayout map;
12976d2a658fSstefano_zampini 
12986d2a658fSstefano_zampini       ierr = MatGetLayouts(T,NULL,&map);CHKERRQ(ierr);
12996d2a658fSstefano_zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
13001c92d2d0Sstefano_zampini       hypre_ParCSRMatrixColStarts(parcsr) = (HYPRE_Int*)(map->range + myid);
13016d2a658fSstefano_zampini     }
13026d2a658fSstefano_zampini     if (!hypre_ParCSRMatrixOwnsRowStarts(parcsr)) {
13036d2a658fSstefano_zampini       PetscLayout map;
13046d2a658fSstefano_zampini 
13056d2a658fSstefano_zampini       ierr = MatGetLayouts(T,&map,NULL);CHKERRQ(ierr);
13066d2a658fSstefano_zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
13071c92d2d0Sstefano_zampini       hypre_ParCSRMatrixRowStarts(parcsr) = (HYPRE_Int*)(map->range + myid);
13086d2a658fSstefano_zampini     }
1309978814f1SStefano Zampini     /* prevent from freeing the pointer */
1310978814f1SStefano Zampini     if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE;
1311225daaf8SStefano Zampini     *A   = T;
1312978814f1SStefano Zampini     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1313978814f1SStefano Zampini     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1314bb4689ddSStefano Zampini   } else if (isaij) {
1315bb4689ddSStefano Zampini     if (copymode != PETSC_OWN_POINTER) {
1316225daaf8SStefano Zampini       /* prevent from freeing the pointer */
1317225daaf8SStefano Zampini       hA->inner_free = PETSC_FALSE;
1318225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
1319225daaf8SStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
1320225daaf8SStefano Zampini     } else { /* AIJ return type with PETSC_OWN_POINTER */
1321225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1322225daaf8SStefano Zampini       *A   = T;
1323225daaf8SStefano Zampini     }
1324bb4689ddSStefano Zampini   } else if (isis) {
13258cfe8d00SStefano Zampini     ierr = MatConvert_HYPRE_IS(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
13268cfe8d00SStefano Zampini     if (copymode != PETSC_OWN_POINTER) hA->inner_free = PETSC_FALSE;
1327bb4689ddSStefano Zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
1328bb4689ddSStefano Zampini   }
1329978814f1SStefano Zampini   PetscFunctionReturn(0);
1330978814f1SStefano Zampini }
1331978814f1SStefano Zampini 
133268ec7858SStefano Zampini static PetscErrorCode MatHYPREGetParCSR_HYPRE(Mat A, hypre_ParCSRMatrix **parcsr)
1333dd9c0a25Sstefano_zampini {
1334dd9c0a25Sstefano_zampini   Mat_HYPRE*            hA = (Mat_HYPRE*)A->data;
1335dd9c0a25Sstefano_zampini   HYPRE_Int             type;
1336dd9c0a25Sstefano_zampini   PetscErrorCode        ierr;
1337dd9c0a25Sstefano_zampini 
1338dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1339dd9c0a25Sstefano_zampini   if (!hA->ij) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"HYPRE_IJMatrix not present");
1340dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
1341dd9c0a25Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"HYPRE_IJMatrix is not of type HYPRE_PARCSR");
1342dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)parcsr));
1343dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1344dd9c0a25Sstefano_zampini }
1345dd9c0a25Sstefano_zampini 
1346dd9c0a25Sstefano_zampini /*
1347dd9c0a25Sstefano_zampini    MatHYPREGetParCSR - Gets the pointer to the ParCSR matrix
1348dd9c0a25Sstefano_zampini 
1349dd9c0a25Sstefano_zampini    Not collective
1350dd9c0a25Sstefano_zampini 
1351dd9c0a25Sstefano_zampini    Input Parameters:
1352dd9c0a25Sstefano_zampini +  A  - the MATHYPRE object
1353dd9c0a25Sstefano_zampini 
1354dd9c0a25Sstefano_zampini    Output Parameter:
1355dd9c0a25Sstefano_zampini .  parcsr  - the pointer to the hypre_ParCSRMatrix
1356dd9c0a25Sstefano_zampini 
1357dd9c0a25Sstefano_zampini    Level: intermediate
1358dd9c0a25Sstefano_zampini 
1359dd9c0a25Sstefano_zampini .seealso: MatHYPRE, PetscCopyMode
1360dd9c0a25Sstefano_zampini */
1361dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR(Mat A, hypre_ParCSRMatrix **parcsr)
1362dd9c0a25Sstefano_zampini {
1363dd9c0a25Sstefano_zampini   PetscErrorCode ierr;
1364dd9c0a25Sstefano_zampini 
1365dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1366dd9c0a25Sstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1367dd9c0a25Sstefano_zampini   PetscValidType(A,1);
1368dd9c0a25Sstefano_zampini   ierr = PetscUseMethod(A,"MatHYPREGetParCSR_C",(Mat,hypre_ParCSRMatrix**),(A,parcsr));CHKERRQ(ierr);
1369dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1370dd9c0a25Sstefano_zampini }
1371dd9c0a25Sstefano_zampini 
137268ec7858SStefano Zampini static PetscErrorCode MatMissingDiagonal_HYPRE(Mat A, PetscBool *missing, PetscInt *dd)
137368ec7858SStefano Zampini {
137468ec7858SStefano Zampini   hypre_ParCSRMatrix *parcsr;
137568ec7858SStefano Zampini   hypre_CSRMatrix    *ha;
137668ec7858SStefano Zampini   PetscInt           rst;
137768ec7858SStefano Zampini   PetscErrorCode     ierr;
137868ec7858SStefano Zampini 
137968ec7858SStefano Zampini   PetscFunctionBegin;
138068ec7858SStefano Zampini   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented with non-square diagonal blocks");
138168ec7858SStefano Zampini   ierr = MatGetOwnershipRange(A,&rst,NULL);CHKERRQ(ierr);
138268ec7858SStefano Zampini   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
138368ec7858SStefano Zampini   if (missing) *missing = PETSC_FALSE;
138468ec7858SStefano Zampini   if (dd) *dd = -1;
138568ec7858SStefano Zampini   ha = hypre_ParCSRMatrixDiag(parcsr);
138668ec7858SStefano Zampini   if (ha) {
138768299464SStefano Zampini     PetscInt  size,i;
138868299464SStefano Zampini     HYPRE_Int *ii,*jj;
138968ec7858SStefano Zampini 
139068ec7858SStefano Zampini     size = hypre_CSRMatrixNumRows(ha);
139168ec7858SStefano Zampini     ii   = hypre_CSRMatrixI(ha);
139268ec7858SStefano Zampini     jj   = hypre_CSRMatrixJ(ha);
139368ec7858SStefano Zampini     for (i = 0; i < size; i++) {
139468ec7858SStefano Zampini       PetscInt  j;
139568ec7858SStefano Zampini       PetscBool found = PETSC_FALSE;
139668ec7858SStefano Zampini 
139768ec7858SStefano Zampini       for (j = ii[i]; j < ii[i+1] && !found; j++)
139868ec7858SStefano Zampini         found = (jj[j] == i) ? PETSC_TRUE : PETSC_FALSE;
139968ec7858SStefano Zampini 
140068ec7858SStefano Zampini       if (!found) {
140168ec7858SStefano Zampini         PetscInfo1(A,"Matrix is missing local diagonal entry %D\n",i);
140268ec7858SStefano Zampini         if (missing) *missing = PETSC_TRUE;
140368ec7858SStefano Zampini         if (dd) *dd = i+rst;
140468ec7858SStefano Zampini         PetscFunctionReturn(0);
140568ec7858SStefano Zampini       }
140668ec7858SStefano Zampini     }
140768ec7858SStefano Zampini     if (!size) {
140868ec7858SStefano Zampini       PetscInfo(A,"Matrix has no diagonal entries therefore is missing diagonal\n");
140968ec7858SStefano Zampini       if (missing) *missing = PETSC_TRUE;
141068ec7858SStefano Zampini       if (dd) *dd = rst;
141168ec7858SStefano Zampini     }
141268ec7858SStefano Zampini   } else {
141368ec7858SStefano Zampini     PetscInfo(A,"Matrix has no diagonal entries therefore is missing diagonal\n");
141468ec7858SStefano Zampini     if (missing) *missing = PETSC_TRUE;
141568ec7858SStefano Zampini     if (dd) *dd = rst;
141668ec7858SStefano Zampini   }
141768ec7858SStefano Zampini   PetscFunctionReturn(0);
141868ec7858SStefano Zampini }
141968ec7858SStefano Zampini 
142068ec7858SStefano Zampini static PetscErrorCode MatScale_HYPRE(Mat A, PetscScalar s)
142168ec7858SStefano Zampini {
142268ec7858SStefano Zampini   hypre_ParCSRMatrix *parcsr;
142368ec7858SStefano Zampini   hypre_CSRMatrix    *ha;
142468ec7858SStefano Zampini   PetscErrorCode     ierr;
142568ec7858SStefano Zampini 
142668ec7858SStefano Zampini   PetscFunctionBegin;
142768ec7858SStefano Zampini   ierr  = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
142868ec7858SStefano Zampini   /* diagonal part */
142968ec7858SStefano Zampini   ha = hypre_ParCSRMatrixDiag(parcsr);
143068ec7858SStefano Zampini   if (ha) {
143168299464SStefano Zampini     PetscInt    size,i;
143268299464SStefano Zampini     HYPRE_Int   *ii;
143368ec7858SStefano Zampini     PetscScalar *a;
143468ec7858SStefano Zampini 
143568ec7858SStefano Zampini     size = hypre_CSRMatrixNumRows(ha);
143668ec7858SStefano Zampini     a    = hypre_CSRMatrixData(ha);
143768ec7858SStefano Zampini     ii   = hypre_CSRMatrixI(ha);
143868ec7858SStefano Zampini     for (i = 0; i < ii[size]; i++) a[i] *= s;
143968ec7858SStefano Zampini   }
144068ec7858SStefano Zampini   /* offdiagonal part */
144168ec7858SStefano Zampini   ha = hypre_ParCSRMatrixOffd(parcsr);
144268ec7858SStefano Zampini   if (ha) {
144368299464SStefano Zampini     PetscInt    size,i;
144468299464SStefano Zampini     HYPRE_Int   *ii;
144568ec7858SStefano Zampini     PetscScalar *a;
144668ec7858SStefano Zampini 
144768ec7858SStefano Zampini     size = hypre_CSRMatrixNumRows(ha);
144868ec7858SStefano Zampini     a    = hypre_CSRMatrixData(ha);
144968ec7858SStefano Zampini     ii   = hypre_CSRMatrixI(ha);
145068ec7858SStefano Zampini     for (i = 0; i < ii[size]; i++) a[i] *= s;
145168ec7858SStefano Zampini   }
145268ec7858SStefano Zampini   PetscFunctionReturn(0);
145368ec7858SStefano Zampini }
145468ec7858SStefano Zampini 
145568ec7858SStefano Zampini static PetscErrorCode MatZeroRowsColumns_HYPRE(Mat A, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
145668ec7858SStefano Zampini {
145768ec7858SStefano Zampini   hypre_ParCSRMatrix *parcsr;
145868299464SStefano Zampini   HYPRE_Int          *lrows;
145968299464SStefano Zampini   PetscInt           rst,ren,i;
146068ec7858SStefano Zampini   PetscErrorCode     ierr;
146168ec7858SStefano Zampini 
146268ec7858SStefano Zampini   PetscFunctionBegin;
146368ec7858SStefano Zampini   if (x || b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");
146468ec7858SStefano Zampini   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
146568ec7858SStefano Zampini   ierr = PetscMalloc1(numRows,&lrows);CHKERRQ(ierr);
146668ec7858SStefano Zampini   ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
146768ec7858SStefano Zampini   for (i=0;i<numRows;i++) {
146868ec7858SStefano Zampini     if (rows[i] < rst || rows[i] >= ren)
146968ec7858SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Non-local rows not yet supported");
147068ec7858SStefano Zampini     lrows[i] = rows[i] - rst;
147168ec7858SStefano Zampini   }
147268ec7858SStefano Zampini   PetscStackCallStandard(hypre_ParCSRMatrixEliminateRowsCols,(parcsr,numRows,lrows));
147368ec7858SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
147468ec7858SStefano Zampini   PetscFunctionReturn(0);
147568ec7858SStefano Zampini }
147668ec7858SStefano Zampini 
1477c69f721fSFande Kong static PetscErrorCode MatZeroEntries_HYPRE_CSRMatrix(hypre_CSRMatrix *ha)
1478c69f721fSFande Kong {
1479c69f721fSFande Kong   PetscErrorCode      ierr;
1480c69f721fSFande Kong 
1481c69f721fSFande Kong   PetscFunctionBegin;
1482c69f721fSFande Kong   if (ha) {
1483c69f721fSFande Kong     HYPRE_Int     *ii, size;
1484c69f721fSFande Kong     HYPRE_Complex *a;
1485c69f721fSFande Kong 
1486c69f721fSFande Kong     size = hypre_CSRMatrixNumRows(ha);
1487c69f721fSFande Kong     a    = hypre_CSRMatrixData(ha);
1488c69f721fSFande Kong     ii   = hypre_CSRMatrixI(ha);
1489c69f721fSFande Kong 
1490c69f721fSFande Kong     if (a) { ierr = PetscMemzero(a,(ii[size])*sizeof(HYPRE_Complex));CHKERRQ(ierr); }
1491c69f721fSFande Kong   }
1492c69f721fSFande Kong   PetscFunctionReturn(0);
1493c69f721fSFande Kong }
1494c69f721fSFande Kong 
1495c69f721fSFande Kong PetscErrorCode MatZeroEntries_HYPRE(Mat A)
1496c69f721fSFande Kong {
1497c69f721fSFande Kong   hypre_ParCSRMatrix  *parcsr;
1498c69f721fSFande Kong   PetscErrorCode      ierr;
1499c69f721fSFande Kong 
1500c69f721fSFande Kong   PetscFunctionBegin;
1501c69f721fSFande Kong   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1502c69f721fSFande Kong   /* diagonal part */
1503c69f721fSFande Kong   ierr = MatZeroEntries_HYPRE_CSRMatrix(hypre_ParCSRMatrixDiag(parcsr));CHKERRQ(ierr);
1504c69f721fSFande Kong   /* off-diagonal part */
1505c69f721fSFande Kong   ierr = MatZeroEntries_HYPRE_CSRMatrix(hypre_ParCSRMatrixOffd(parcsr));CHKERRQ(ierr);
1506c69f721fSFande Kong   PetscFunctionReturn(0);
1507c69f721fSFande Kong }
1508c69f721fSFande Kong 
1509c69f721fSFande Kong 
1510c69f721fSFande Kong static PetscErrorCode MatZeroRows_HYPRE_CSRMatrix(hypre_CSRMatrix *hA,PetscInt N,const PetscInt rows[],PetscScalar diag)
1511c69f721fSFande Kong {
1512c69f721fSFande Kong   PetscInt        ii, jj, ibeg, iend, irow;
1513c69f721fSFande Kong   PetscInt        *i, *j;
1514c69f721fSFande Kong   PetscScalar     *a;
1515c69f721fSFande Kong 
1516c69f721fSFande Kong   PetscFunctionBegin;
1517c69f721fSFande Kong 
1518c69f721fSFande Kong   if (!hA) PetscFunctionReturn(0);
1519c69f721fSFande Kong 
1520*08defe43SFande Kong   i = (PetscInt*) hypre_CSRMatrixI(hA);
1521*08defe43SFande Kong   j = (PetscInt*) hypre_CSRMatrixJ(hA);
1522c69f721fSFande Kong   a = hypre_CSRMatrixData(hA);
1523c69f721fSFande Kong 
1524c69f721fSFande Kong   for (ii = 0; ii < N; ii++) {
1525c69f721fSFande Kong     irow = rows[ii];
1526c69f721fSFande Kong     ibeg = i[irow];
1527c69f721fSFande Kong     iend = i[irow+1];
1528c69f721fSFande Kong     for (jj = ibeg; jj < iend; jj++)
1529c69f721fSFande Kong       if (j[jj] == irow) a[jj] = diag;
1530c69f721fSFande Kong       else a[jj] = 0.0;
1531c69f721fSFande Kong    }
1532c69f721fSFande Kong 
1533c69f721fSFande Kong    PetscFunctionReturn(0);
1534c69f721fSFande Kong }
1535c69f721fSFande Kong 
1536c69f721fSFande Kong PetscErrorCode MatZeroRows_HYPRE(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1537c69f721fSFande Kong {
1538c69f721fSFande Kong   hypre_ParCSRMatrix  *parcsr;
1539c69f721fSFande Kong   PetscInt            *lrows,len;
1540c69f721fSFande Kong   PetscErrorCode      ierr;
1541c69f721fSFande Kong 
1542c69f721fSFande Kong   PetscFunctionBegin;
1543c69f721fSFande Kong   if (x || b) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Does not support to modify the solution and the right hand size \n");
1544c69f721fSFande Kong   /* retrieve the internal matrix */
1545c69f721fSFande Kong   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1546c69f721fSFande Kong   /* get locally owned rows */
1547c69f721fSFande Kong   ierr = MatZeroRowsMapLocal_Private(A,N,rows,&len,&lrows);CHKERRQ(ierr);
1548c69f721fSFande Kong   /* zero diagonal part */
1549c69f721fSFande Kong   ierr = MatZeroRows_HYPRE_CSRMatrix(hypre_ParCSRMatrixDiag(parcsr),len,lrows,diag);CHKERRQ(ierr);
1550c69f721fSFande Kong   /* zero off-diagonal part */
1551c69f721fSFande Kong   ierr = MatZeroRows_HYPRE_CSRMatrix(hypre_ParCSRMatrixOffd(parcsr),len,lrows,0.0);CHKERRQ(ierr);
1552c69f721fSFande Kong 
1553c69f721fSFande Kong   ierr = PetscFree(lrows);CHKERRQ(ierr);
1554c69f721fSFande Kong   PetscFunctionReturn(0);
1555c69f721fSFande Kong }
1556c69f721fSFande Kong 
1557c69f721fSFande Kong 
1558c69f721fSFande Kong PetscErrorCode MatAssemblyBegin_HYPRE(Mat mat,MatAssemblyType mode)
1559c69f721fSFande Kong {
1560c69f721fSFande Kong   PetscErrorCode ierr;
1561c69f721fSFande Kong 
1562c69f721fSFande Kong   PetscFunctionBegin;
1563c69f721fSFande Kong   if (mat->nooffprocentries) PetscFunctionReturn(0);
1564c69f721fSFande Kong 
1565c69f721fSFande Kong   ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr);
1566c69f721fSFande Kong   PetscFunctionReturn(0);
1567c69f721fSFande Kong }
1568c69f721fSFande Kong 
1569c69f721fSFande Kong PetscErrorCode MatGetRow_HYPRE(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1570c69f721fSFande Kong {
1571c69f721fSFande Kong   hypre_ParCSRMatrix  *parcsr;
1572c69f721fSFande Kong   PetscErrorCode      ierr;
1573c69f721fSFande Kong 
1574c69f721fSFande Kong   PetscFunctionBegin;
1575c69f721fSFande Kong   /* retrieve the internal matrix */
1576c69f721fSFande Kong   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1577c69f721fSFande Kong   /* call HYPRE API */
1578*08defe43SFande Kong   PetscStackCallStandard(HYPRE_ParCSRMatrixGetRow,(parcsr,row,(HYPRE_Int*)nz,(HYPRE_Int**)idx,v));
1579c69f721fSFande Kong   PetscFunctionReturn(0);
1580c69f721fSFande Kong }
1581c69f721fSFande Kong 
1582c69f721fSFande Kong PetscErrorCode MatRestoreRow_HYPRE(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1583c69f721fSFande Kong {
1584c69f721fSFande Kong   hypre_ParCSRMatrix  *parcsr;
1585c69f721fSFande Kong   PetscErrorCode      ierr;
1586c69f721fSFande Kong 
1587c69f721fSFande Kong   PetscFunctionBegin;
1588c69f721fSFande Kong   /* retrieve the internal matrix */
1589c69f721fSFande Kong   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1590c69f721fSFande Kong   /* call HYPRE API */
1591*08defe43SFande Kong   PetscStackCallStandard(HYPRE_ParCSRMatrixRestoreRow,(parcsr,row,(HYPRE_Int*)nz,(HYPRE_Int**)idx,v));
1592c69f721fSFande Kong   PetscFunctionReturn(0);
1593c69f721fSFande Kong }
1594c69f721fSFande Kong 
1595c69f721fSFande Kong 
1596c69f721fSFande Kong 
1597c69f721fSFande Kong PetscErrorCode MatGetValues_HYPRE(Mat A,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
1598c69f721fSFande Kong {
1599c69f721fSFande Kong   HYPRE_IJMatrix     *hIJ = (HYPRE_IJMatrix*)A->data;
1600c69f721fSFande Kong   PetscErrorCode      ierr;
1601c69f721fSFande Kong   PetscInt            i;
1602c69f721fSFande Kong   PetscFunctionBegin;
1603c69f721fSFande Kong   if (!m || !n) PetscFunctionReturn(0);
1604c69f721fSFande Kong 
1605c69f721fSFande Kong   /* Ignore negative row indices
1606c69f721fSFande Kong    * And negative column indices should be automatically ignored in hypre
1607c69f721fSFande Kong    * */
1608c69f721fSFande Kong   for (i=0; i<m; i++)
1609*08defe43SFande Kong     if (idxm[i] >= 0) PetscStackCallStandard(HYPRE_IJMatrixGetValues,(*hIJ,1,(HYPRE_Int*)&n,(HYPRE_Int*)&idxm[i],(HYPRE_Int*)idxn,&v[i*n]));
1610c69f721fSFande Kong 
1611c69f721fSFande Kong   PetscFunctionReturn(0);
1612c69f721fSFande Kong }
1613c69f721fSFande Kong 
1614c69f721fSFande Kong 
1615a055b5aaSBarry Smith /*MC
1616a055b5aaSBarry Smith    MATHYPRE - MATHYPRE = "hypre" - A matrix type to be used for sequential and parallel sparse matrices
1617a055b5aaSBarry Smith           based on the hypre IJ interface.
1618a055b5aaSBarry Smith 
1619a055b5aaSBarry Smith    Level: intermediate
1620a055b5aaSBarry Smith 
1621a055b5aaSBarry Smith .seealso: MatCreate()
1622a055b5aaSBarry Smith M*/
1623a055b5aaSBarry Smith 
162463c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
162563c07aadSStefano Zampini {
162663c07aadSStefano Zampini   Mat_HYPRE      *hB;
162763c07aadSStefano Zampini   PetscErrorCode ierr;
162863c07aadSStefano Zampini 
162963c07aadSStefano Zampini   PetscFunctionBegin;
163063c07aadSStefano Zampini   ierr           = PetscNewLog(B,&hB);CHKERRQ(ierr);
1631978814f1SStefano Zampini   hB->inner_free = PETSC_TRUE;
1632c69f721fSFande Kong   hB->available  = PETSC_TRUE;
1633c69f721fSFande Kong   hB->size       = 0;
1634c69f721fSFande Kong   hB->array      = NULL;
1635978814f1SStefano Zampini 
163663c07aadSStefano Zampini   B->data       = (void*)hB;
163763c07aadSStefano Zampini   B->rmap->bs   = 1;
163863c07aadSStefano Zampini   B->assembled  = PETSC_FALSE;
163963c07aadSStefano Zampini 
1640de393100SStefano Zampini   ierr = PetscMemzero(B->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
164163c07aadSStefano Zampini   B->ops->mult            = MatMult_HYPRE;
164263c07aadSStefano Zampini   B->ops->multtranspose   = MatMultTranspose_HYPRE;
164363c07aadSStefano Zampini   B->ops->setup           = MatSetUp_HYPRE;
164463c07aadSStefano Zampini   B->ops->destroy         = MatDestroy_HYPRE;
1645c69f721fSFande Kong 
1646c69f721fSFande Kong   /* build cache for off array entries formed */
1647c69f721fSFande Kong   ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);CHKERRQ(ierr);
164863c07aadSStefano Zampini   B->ops->assemblyend     = MatAssemblyEnd_HYPRE;
1649c69f721fSFande Kong   B->ops->assemblybegin   = MatAssemblyBegin_HYPRE;
1650c69f721fSFande Kong 
1651cd8bc7baSStefano Zampini   B->ops->ptap            = MatPtAP_HYPRE_HYPRE;
1652d501dc42Sstefano_zampini   B->ops->matmult         = MatMatMult_HYPRE_HYPRE;
1653d975228cSstefano_zampini   B->ops->setvalues       = MatSetValues_HYPRE;
165468ec7858SStefano Zampini   B->ops->missingdiagonal = MatMissingDiagonal_HYPRE;
165568ec7858SStefano Zampini   B->ops->scale           = MatScale_HYPRE;
165668ec7858SStefano Zampini   B->ops->zerorowscolumns = MatZeroRowsColumns_HYPRE;
1657c69f721fSFande Kong   B->ops->zeroentries     = MatZeroEntries_HYPRE;
1658c69f721fSFande Kong   B->ops->zerorows        = MatZeroRows_HYPRE;
1659c69f721fSFande Kong   B->ops->getrow          = MatGetRow_HYPRE;
1660c69f721fSFande Kong   B->ops->restorerow      = MatRestoreRow_HYPRE;
1661c69f721fSFande Kong   B->ops->getvalues       = MatGetValues_HYPRE;
166263c07aadSStefano Zampini 
166363c07aadSStefano Zampini   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
166463c07aadSStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
166563c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
16662df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr);
1667c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1668c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1669d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPRESetPreallocation_C",MatHYPRESetPreallocation_HYPRE);CHKERRQ(ierr);
1670dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPREGetParCSR_C",MatHYPREGetParCSR_HYPRE);CHKERRQ(ierr);
167163c07aadSStefano Zampini   PetscFunctionReturn(0);
167263c07aadSStefano Zampini }
167363c07aadSStefano Zampini 
1674225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void *ptr)
1675225daaf8SStefano Zampini {
1676225daaf8SStefano Zampini    PetscFunctionBegin;
1677e6de0934SSatish Balay    hypre_TFree(ptr,HYPRE_MEMORY_HOST);
1678225daaf8SStefano Zampini    PetscFunctionReturn(0);
1679225daaf8SStefano Zampini }
1680