xref: /petsc/src/mat/impls/hypre/mhypre.c (revision 75d48cdbf3d989106dd25a7997baa36b63209129)
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 
16*75d48cdbSStefano Zampini PETSC_INTERN PetscErrorCode MatPtAP_IS_XAIJ(Mat,Mat,MatReuse,PetscReal,Mat*);
17*75d48cdbSStefano Zampini 
1863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*);
1963c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix);
2063c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
2163c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
2263c07aadSStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool);
23225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void*);
2463c07aadSStefano Zampini 
2563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij)
2663c07aadSStefano Zampini {
2763c07aadSStefano Zampini   PetscErrorCode ierr;
2863c07aadSStefano Zampini   PetscInt       i,n_d,n_o;
2963c07aadSStefano Zampini   const PetscInt *ia_d,*ia_o;
3063c07aadSStefano Zampini   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
3163c07aadSStefano Zampini   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
3263c07aadSStefano Zampini 
3363c07aadSStefano Zampini   PetscFunctionBegin;
3463c07aadSStefano Zampini   if (A_d) { /* determine number of nonzero entries in local diagonal part */
3563c07aadSStefano Zampini     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
3663c07aadSStefano Zampini     if (done_d) {
3763c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr);
3863c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
3963c07aadSStefano Zampini         nnz_d[i] = ia_d[i+1] - ia_d[i];
4063c07aadSStefano Zampini       }
4163c07aadSStefano Zampini     }
4263c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
4363c07aadSStefano Zampini   }
4463c07aadSStefano Zampini   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
4563c07aadSStefano Zampini     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
4663c07aadSStefano Zampini     if (done_o) {
4763c07aadSStefano Zampini       ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr);
4863c07aadSStefano Zampini       for (i=0; i<n_o; i++) {
4963c07aadSStefano Zampini         nnz_o[i] = ia_o[i+1] - ia_o[i];
5063c07aadSStefano Zampini       }
5163c07aadSStefano Zampini     }
5263c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
5363c07aadSStefano Zampini   }
5463c07aadSStefano Zampini   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
5563c07aadSStefano Zampini     if (!done_o) { /* only diagonal part */
5663c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
5763c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
5863c07aadSStefano Zampini         nnz_o[i] = 0;
5963c07aadSStefano Zampini       }
6063c07aadSStefano Zampini     }
6163c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o));
6263c07aadSStefano Zampini     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
6363c07aadSStefano Zampini     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
6463c07aadSStefano Zampini   }
6563c07aadSStefano Zampini   PetscFunctionReturn(0);
6663c07aadSStefano Zampini }
6763c07aadSStefano Zampini 
6863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA)
6963c07aadSStefano Zampini {
7063c07aadSStefano Zampini   PetscErrorCode ierr;
7163c07aadSStefano Zampini   PetscInt       rstart,rend,cstart,cend;
7263c07aadSStefano Zampini 
7363c07aadSStefano Zampini   PetscFunctionBegin;
74d975228cSstefano_zampini   ierr   = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
75d975228cSstefano_zampini   ierr   = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
7663c07aadSStefano Zampini   rstart = A->rmap->rstart;
7763c07aadSStefano Zampini   rend   = A->rmap->rend;
7863c07aadSStefano Zampini   cstart = A->cmap->rstart;
7963c07aadSStefano Zampini   cend   = A->cmap->rend;
8063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
8163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
8263c07aadSStefano Zampini   {
8363c07aadSStefano Zampini     PetscBool      same;
8463c07aadSStefano Zampini     Mat            A_d,A_o;
8563c07aadSStefano Zampini     const PetscInt *colmap;
8663c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
8763c07aadSStefano Zampini     if (same) {
8863c07aadSStefano Zampini       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
8963c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
9063c07aadSStefano Zampini       PetscFunctionReturn(0);
9163c07aadSStefano Zampini     }
9263c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
9363c07aadSStefano Zampini     if (same) {
9463c07aadSStefano Zampini       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
9563c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
9663c07aadSStefano Zampini       PetscFunctionReturn(0);
9763c07aadSStefano Zampini     }
9863c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
9963c07aadSStefano Zampini     if (same) {
10063c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
10163c07aadSStefano Zampini       PetscFunctionReturn(0);
10263c07aadSStefano Zampini     }
10363c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
10463c07aadSStefano Zampini     if (same) {
10563c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
10663c07aadSStefano Zampini       PetscFunctionReturn(0);
10763c07aadSStefano Zampini     }
10863c07aadSStefano Zampini   }
10963c07aadSStefano Zampini   PetscFunctionReturn(0);
11063c07aadSStefano Zampini }
11163c07aadSStefano Zampini 
11263c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij)
11363c07aadSStefano Zampini {
11463c07aadSStefano Zampini   PetscErrorCode    ierr;
11563c07aadSStefano Zampini   PetscInt          i,rstart,rend,ncols,nr,nc;
11663c07aadSStefano Zampini   const PetscScalar *values;
11763c07aadSStefano Zampini   const PetscInt    *cols;
11863c07aadSStefano Zampini   PetscBool         flg;
11963c07aadSStefano Zampini 
12063c07aadSStefano Zampini   PetscFunctionBegin;
12163c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
12263c07aadSStefano Zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
12363c07aadSStefano Zampini   if (flg && nr == nc) {
12463c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
12563c07aadSStefano Zampini     PetscFunctionReturn(0);
12663c07aadSStefano Zampini   }
12763c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
12863c07aadSStefano Zampini   if (flg) {
12963c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
13063c07aadSStefano Zampini     PetscFunctionReturn(0);
13163c07aadSStefano Zampini   }
13263c07aadSStefano Zampini 
13363c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
13463c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
13563c07aadSStefano Zampini   for (i=rstart; i<rend; i++) {
13663c07aadSStefano Zampini     ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
137e3977e59Sstefano_zampini     if (ncols) {
13863c07aadSStefano Zampini       PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
139e3977e59Sstefano_zampini     }
14063c07aadSStefano Zampini     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
14163c07aadSStefano Zampini   }
14263c07aadSStefano Zampini   PetscFunctionReturn(0);
14363c07aadSStefano Zampini }
14463c07aadSStefano Zampini 
14563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij)
14663c07aadSStefano Zampini {
14763c07aadSStefano Zampini   PetscErrorCode        ierr;
14863c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag = (Mat_SeqAIJ*)A->data;
14958968eb6SStefano Zampini   HYPRE_Int             type;
15063c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
15163c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
15263c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag;
15363c07aadSStefano Zampini 
15463c07aadSStefano Zampini   PetscFunctionBegin;
15563c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
156ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
157ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
158ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
15963c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
16063c07aadSStefano Zampini   /*
16163c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
16263c07aadSStefano Zampini   */
16363c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));
16463c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));
16563c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
166ea9daf28SStefano Zampini 
167ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
16863c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
16963c07aadSStefano Zampini   PetscFunctionReturn(0);
17063c07aadSStefano Zampini }
17163c07aadSStefano Zampini 
17263c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
17363c07aadSStefano Zampini {
17463c07aadSStefano Zampini   PetscErrorCode        ierr;
17563c07aadSStefano Zampini   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
17663c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag,*poffd;
17763c07aadSStefano Zampini   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
17858968eb6SStefano Zampini   HYPRE_Int             type;
17963c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
18063c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
18163c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag,*hoffd;
18263c07aadSStefano Zampini 
18363c07aadSStefano Zampini   PetscFunctionBegin;
18463c07aadSStefano Zampini   pdiag = (Mat_SeqAIJ*) pA->A->data;
18563c07aadSStefano Zampini   poffd = (Mat_SeqAIJ*) pA->B->data;
18663c07aadSStefano Zampini   /* cstart is only valid for square MPIAIJ layed out in the usual way */
18763c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
18863c07aadSStefano Zampini 
18963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
190ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
191ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
192ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
19363c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
19463c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
19563c07aadSStefano Zampini 
19663c07aadSStefano Zampini   /*
19763c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
19863c07aadSStefano Zampini   */
19963c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
20063c07aadSStefano Zampini   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
20163c07aadSStefano Zampini   jj  = (PetscInt*)hdiag->j;
20263c07aadSStefano Zampini   pjj = (PetscInt*)pdiag->j;
20363c07aadSStefano Zampini   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
20463c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
20563c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
20663c07aadSStefano Zampini   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
20763c07aadSStefano Zampini      If we hacked a hypre a bit more we might be able to avoid this step */
20863c07aadSStefano Zampini   jj  = (PetscInt*) hoffd->j;
20963c07aadSStefano Zampini   pjj = (PetscInt*) poffd->j;
21063c07aadSStefano Zampini   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
21163c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
21263c07aadSStefano Zampini 
213ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
21463c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
21563c07aadSStefano Zampini   PetscFunctionReturn(0);
21663c07aadSStefano Zampini }
21763c07aadSStefano Zampini 
2182df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B)
2192df22349SStefano Zampini {
2202df22349SStefano Zampini   Mat_HYPRE*             mhA = (Mat_HYPRE*)(A->data);
2212df22349SStefano Zampini   Mat                    lA;
2222df22349SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
2232df22349SStefano Zampini   IS                     is;
2242df22349SStefano Zampini   hypre_ParCSRMatrix     *hA;
2252df22349SStefano Zampini   hypre_CSRMatrix        *hdiag,*hoffd;
2262df22349SStefano Zampini   MPI_Comm               comm;
2272df22349SStefano Zampini   PetscScalar            *hdd,*hod,*aa,*data;
2287d968826Sstefano_zampini   HYPRE_Int              *col_map_offd,*hdi,*hdj,*hoi,*hoj;
2292df22349SStefano Zampini   PetscInt               *ii,*jj,*iptr,*jptr;
2302df22349SStefano Zampini   PetscInt               cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N;
23158968eb6SStefano Zampini   HYPRE_Int              type;
2322df22349SStefano Zampini   PetscErrorCode         ierr;
2332df22349SStefano Zampini 
2342df22349SStefano Zampini   PetscFunctionBegin;
235a1787963SStefano Zampini   comm = PetscObjectComm((PetscObject)A);
2362df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type));
2372df22349SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
2382df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA));
2392df22349SStefano Zampini   M     = hypre_ParCSRMatrixGlobalNumRows(hA);
2402df22349SStefano Zampini   N     = hypre_ParCSRMatrixGlobalNumCols(hA);
2412df22349SStefano Zampini   str   = hypre_ParCSRMatrixFirstRowIndex(hA);
2422df22349SStefano Zampini   stc   = hypre_ParCSRMatrixFirstColDiag(hA);
2432df22349SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(hA);
2442df22349SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(hA);
2452df22349SStefano Zampini   dr    = hypre_CSRMatrixNumRows(hdiag);
2462df22349SStefano Zampini   dc    = hypre_CSRMatrixNumCols(hdiag);
2472df22349SStefano Zampini   nnz   = hypre_CSRMatrixNumNonzeros(hdiag);
2482df22349SStefano Zampini   hdi   = hypre_CSRMatrixI(hdiag);
2492df22349SStefano Zampini   hdj   = hypre_CSRMatrixJ(hdiag);
2502df22349SStefano Zampini   hdd   = hypre_CSRMatrixData(hdiag);
2512df22349SStefano Zampini   oc    = hypre_CSRMatrixNumCols(hoffd);
2522df22349SStefano Zampini   nnz  += hypre_CSRMatrixNumNonzeros(hoffd);
2532df22349SStefano Zampini   hoi   = hypre_CSRMatrixI(hoffd);
2542df22349SStefano Zampini   hoj   = hypre_CSRMatrixJ(hoffd);
2552df22349SStefano Zampini   hod   = hypre_CSRMatrixData(hoffd);
2562df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
2572df22349SStefano Zampini     PetscInt *aux;
2582df22349SStefano Zampini 
2592df22349SStefano Zampini     /* generate l2g maps for rows and cols */
2602df22349SStefano Zampini     ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr);
2612df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
2622df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2632df22349SStefano Zampini     col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
2642df22349SStefano Zampini     ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr);
2652df22349SStefano Zampini     for (i=0; i<dc; i++) aux[i] = i+stc;
2662df22349SStefano Zampini     for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i];
2672df22349SStefano Zampini     ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
2682df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
2692df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2702df22349SStefano Zampini     /* create MATIS object */
2712df22349SStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
2722df22349SStefano Zampini     ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr);
2732df22349SStefano Zampini     ierr = MatSetType(*B,MATIS);CHKERRQ(ierr);
2742df22349SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr);
2752df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
2762df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
2772df22349SStefano Zampini 
2782df22349SStefano Zampini     /* allocate CSR for local matrix */
2792df22349SStefano Zampini     ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr);
2802df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr);
2812df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr);
2822df22349SStefano Zampini   } else {
2832df22349SStefano Zampini     PetscInt  nr;
2842df22349SStefano Zampini     PetscBool done;
2852df22349SStefano Zampini     ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr);
2862df22349SStefano Zampini     ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr);
2872df22349SStefano 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);
2882df22349SStefano 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);
2892df22349SStefano Zampini     ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr);
2902df22349SStefano Zampini   }
2912df22349SStefano Zampini   /* merge local matrices */
2922df22349SStefano Zampini   ii   = iptr;
2932df22349SStefano Zampini   jj   = jptr;
2942df22349SStefano Zampini   aa   = data;
2952df22349SStefano Zampini   *ii  = *(hdi++) + *(hoi++);
2962df22349SStefano Zampini   for (jd=0,jo=0,cum=0; *ii<nnz; cum++) {
2972df22349SStefano Zampini     PetscScalar *aold = aa;
2982df22349SStefano Zampini     PetscInt    *jold = jj,nc = jd+jo;
2992df22349SStefano Zampini     for (; jd<*hdi; jd++) { *jj++ = *hdj++;      *aa++ = *hdd++; }
3002df22349SStefano Zampini     for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
3012df22349SStefano Zampini     *(++ii) = *(hdi++) + *(hoi++);
3022df22349SStefano Zampini     ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr);
3032df22349SStefano Zampini   }
3042df22349SStefano Zampini   for (; cum<dr; cum++) *(++ii) = nnz;
3052df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
306a033916dSStefano Zampini     Mat_SeqAIJ* a;
307a033916dSStefano Zampini 
3082df22349SStefano Zampini     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr);
3092df22349SStefano Zampini     ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr);
310a033916dSStefano Zampini     /* hack SeqAIJ */
311a033916dSStefano Zampini     a          = (Mat_SeqAIJ*)(lA->data);
312a033916dSStefano Zampini     a->free_a  = PETSC_TRUE;
313a033916dSStefano Zampini     a->free_ij = PETSC_TRUE;
3142df22349SStefano Zampini     ierr = MatDestroy(&lA);CHKERRQ(ierr);
3152df22349SStefano Zampini   }
3162df22349SStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3172df22349SStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3182df22349SStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
3192df22349SStefano Zampini     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
3202df22349SStefano Zampini   }
3212df22349SStefano Zampini   PetscFunctionReturn(0);
3222df22349SStefano Zampini }
3232df22349SStefano Zampini 
32463c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
32563c07aadSStefano Zampini {
32663c07aadSStefano Zampini   Mat_HYPRE      *hB;
32763c07aadSStefano Zampini   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
32863c07aadSStefano Zampini   PetscErrorCode ierr;
32963c07aadSStefano Zampini 
33063c07aadSStefano Zampini   PetscFunctionBegin;
33163c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
33263c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
33363c07aadSStefano Zampini     /* always destroy the old matrix and create a new memory;
33463c07aadSStefano Zampini        hope this does not churn the memory too much. The problem
33563c07aadSStefano Zampini        is I do not know if it is possible to put the matrix back to
33663c07aadSStefano Zampini        its initial state so that we can directly copy the values
33763c07aadSStefano Zampini        the second time through. */
33863c07aadSStefano Zampini     hB = (Mat_HYPRE*)((*B)->data);
33963c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
34063c07aadSStefano Zampini   } else {
34163c07aadSStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
34263c07aadSStefano Zampini     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
34363c07aadSStefano Zampini     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
34463c07aadSStefano Zampini     hB   = (Mat_HYPRE*)((*B)->data);
34563c07aadSStefano Zampini   }
34663c07aadSStefano Zampini   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
34763c07aadSStefano Zampini   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
3484ec6421dSstefano_zampini   (*B)->preallocated = PETSC_TRUE;
34963c07aadSStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
35063c07aadSStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
35163c07aadSStefano Zampini   PetscFunctionReturn(0);
35263c07aadSStefano Zampini }
35363c07aadSStefano Zampini 
354ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
35563c07aadSStefano Zampini {
35663c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
35763c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
35863c07aadSStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
35963c07aadSStefano Zampini   MPI_Comm           comm;
36063c07aadSStefano Zampini   PetscScalar        *da,*oa,*aptr;
36163c07aadSStefano Zampini   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
36263c07aadSStefano Zampini   PetscInt           i,dnnz,onnz,m,n;
36358968eb6SStefano Zampini   HYPRE_Int          type;
36463c07aadSStefano Zampini   PetscMPIInt        size;
36563c07aadSStefano Zampini   PetscErrorCode     ierr;
36663c07aadSStefano Zampini 
36763c07aadSStefano Zampini   PetscFunctionBegin;
36863c07aadSStefano Zampini   comm = PetscObjectComm((PetscObject)A);
36963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
37063c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
37163c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
37263c07aadSStefano Zampini     PetscBool ismpiaij,isseqaij;
37363c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
3744099cc6bSBarry Smith     ierr = PetscObjectBaseTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
37563c07aadSStefano Zampini     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
37663c07aadSStefano Zampini   }
37763c07aadSStefano Zampini   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
37863c07aadSStefano Zampini 
37963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
38063c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(parcsr);
38163c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(parcsr);
38263c07aadSStefano Zampini   m     = hypre_CSRMatrixNumRows(hdiag);
38363c07aadSStefano Zampini   n     = hypre_CSRMatrixNumCols(hdiag);
38463c07aadSStefano Zampini   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
38563c07aadSStefano Zampini   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
386225daaf8SStefano Zampini   if (reuse == MAT_INITIAL_MATRIX) {
38763c07aadSStefano Zampini     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
38863c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
38963c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
390225daaf8SStefano Zampini   } else if (reuse == MAT_REUSE_MATRIX) {
39163c07aadSStefano Zampini     PetscInt  nr;
39263c07aadSStefano Zampini     PetscBool done;
39363c07aadSStefano Zampini     if (size > 1) {
39463c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
39563c07aadSStefano Zampini 
39663c07aadSStefano Zampini       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
39763c07aadSStefano 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);
39863c07aadSStefano 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);
39963c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
40063c07aadSStefano Zampini     } else {
40163c07aadSStefano Zampini       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
40263c07aadSStefano Zampini       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
40363c07aadSStefano 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);
40463c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
40563c07aadSStefano Zampini     }
406225daaf8SStefano Zampini   } else { /* MAT_INPLACE_MATRIX */
4077d968826Sstefano_zampini     dii = (PetscInt*)hypre_CSRMatrixI(hdiag);
4087d968826Sstefano_zampini     djj = (PetscInt*)hypre_CSRMatrixJ(hdiag);
409225daaf8SStefano Zampini     da  = hypre_CSRMatrixData(hdiag);
41063c07aadSStefano Zampini   }
41163c07aadSStefano Zampini   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
41263c07aadSStefano Zampini   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
41363c07aadSStefano Zampini   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
41463c07aadSStefano Zampini   iptr = djj;
41563c07aadSStefano Zampini   aptr = da;
41663c07aadSStefano Zampini   for (i=0; i<m; i++) {
41763c07aadSStefano Zampini     PetscInt nc = dii[i+1]-dii[i];
41863c07aadSStefano Zampini     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
41963c07aadSStefano Zampini     iptr += nc;
42063c07aadSStefano Zampini     aptr += nc;
42163c07aadSStefano Zampini   }
42263c07aadSStefano Zampini   if (size > 1) {
4237d968826Sstefano_zampini     HYPRE_Int *offdj,*coffd;
42463c07aadSStefano Zampini 
425225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
42663c07aadSStefano Zampini       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
42763c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
42863c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
429225daaf8SStefano Zampini     } else if (reuse == MAT_REUSE_MATRIX) {
43063c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
43163c07aadSStefano Zampini       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
43263c07aadSStefano Zampini       PetscBool  done;
43363c07aadSStefano Zampini 
43463c07aadSStefano Zampini       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
43563c07aadSStefano 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);
43663c07aadSStefano 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);
43763c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
438225daaf8SStefano Zampini     } else { /* MAT_INPLACE_MATRIX */
4397d968826Sstefano_zampini       oii = (PetscInt*)hypre_CSRMatrixI(hoffd);
4407d968826Sstefano_zampini       ojj = (PetscInt*)hypre_CSRMatrixJ(hoffd);
441225daaf8SStefano Zampini       oa  = hypre_CSRMatrixData(hoffd);
44263c07aadSStefano Zampini     }
44363c07aadSStefano Zampini     ierr  = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
44463c07aadSStefano Zampini     offdj = hypre_CSRMatrixJ(hoffd);
44563c07aadSStefano Zampini     coffd = hypre_ParCSRMatrixColMapOffd(parcsr);
44663c07aadSStefano Zampini     for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]];
44763c07aadSStefano Zampini     ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr);
44863c07aadSStefano Zampini     iptr = ojj;
44963c07aadSStefano Zampini     aptr = oa;
45063c07aadSStefano Zampini     for (i=0; i<m; i++) {
45163c07aadSStefano Zampini        PetscInt nc = oii[i+1]-oii[i];
45263c07aadSStefano Zampini        ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
45363c07aadSStefano Zampini        iptr += nc;
45463c07aadSStefano Zampini        aptr += nc;
45563c07aadSStefano Zampini     }
456225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
45763c07aadSStefano Zampini       Mat_MPIAIJ *b;
45863c07aadSStefano Zampini       Mat_SeqAIJ *d,*o;
459225daaf8SStefano Zampini 
46041073d71Sstefano_zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,B);CHKERRQ(ierr);
46163c07aadSStefano Zampini       /* hack MPIAIJ */
46263c07aadSStefano Zampini       b          = (Mat_MPIAIJ*)((*B)->data);
46363c07aadSStefano Zampini       d          = (Mat_SeqAIJ*)b->A->data;
46463c07aadSStefano Zampini       o          = (Mat_SeqAIJ*)b->B->data;
46563c07aadSStefano Zampini       d->free_a  = PETSC_TRUE;
46663c07aadSStefano Zampini       d->free_ij = PETSC_TRUE;
46763c07aadSStefano Zampini       o->free_a  = PETSC_TRUE;
46863c07aadSStefano Zampini       o->free_ij = PETSC_TRUE;
469225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
470225daaf8SStefano Zampini       Mat T;
47141073d71Sstefano_zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,&T);CHKERRQ(ierr);
472225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
473225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
474225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
475225daaf8SStefano Zampini       hypre_CSRMatrixI(hoffd)    = NULL;
476225daaf8SStefano Zampini       hypre_CSRMatrixJ(hoffd)    = NULL;
477225daaf8SStefano Zampini       hypre_CSRMatrixData(hoffd) = NULL;
478225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
47963c07aadSStefano Zampini     }
480225daaf8SStefano Zampini   } else {
481225daaf8SStefano Zampini     oii  = NULL;
482225daaf8SStefano Zampini     ojj  = NULL;
483225daaf8SStefano Zampini     oa   = NULL;
484225daaf8SStefano Zampini     if (reuse == MAT_INITIAL_MATRIX) {
48563c07aadSStefano Zampini       Mat_SeqAIJ* b;
48663c07aadSStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
48763c07aadSStefano Zampini       /* hack SeqAIJ */
48863c07aadSStefano Zampini       b          = (Mat_SeqAIJ*)((*B)->data);
48963c07aadSStefano Zampini       b->free_a  = PETSC_TRUE;
49063c07aadSStefano Zampini       b->free_ij = PETSC_TRUE;
491225daaf8SStefano Zampini     } else if (reuse == MAT_INPLACE_MATRIX) {
492225daaf8SStefano Zampini       Mat T;
493225daaf8SStefano Zampini       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,&T);CHKERRQ(ierr);
494225daaf8SStefano Zampini       hypre_CSRMatrixI(hdiag)    = NULL;
495225daaf8SStefano Zampini       hypre_CSRMatrixJ(hdiag)    = NULL;
496225daaf8SStefano Zampini       hypre_CSRMatrixData(hdiag) = NULL;
497225daaf8SStefano Zampini       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
49863c07aadSStefano Zampini     }
499225daaf8SStefano Zampini   }
500225daaf8SStefano Zampini 
501225daaf8SStefano Zampini   /* we have to use hypre_Tfree to free the arrays */
50263c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
503225daaf8SStefano Zampini     void *ptrs[6] = {dii,djj,da,oii,ojj,oa};
504225daaf8SStefano Zampini     const char *names[6] = {"_hypre_csr_dii",
505225daaf8SStefano Zampini                             "_hypre_csr_djj",
506225daaf8SStefano Zampini                             "_hypre_csr_da",
507225daaf8SStefano Zampini                             "_hypre_csr_oii",
508225daaf8SStefano Zampini                             "_hypre_csr_ojj",
509225daaf8SStefano Zampini                             "_hypre_csr_oa"};
510225daaf8SStefano Zampini     for (i=0; i<6; i++) {
511225daaf8SStefano Zampini       PetscContainer c;
512225daaf8SStefano Zampini 
513225daaf8SStefano Zampini       ierr = PetscContainerCreate(comm,&c);CHKERRQ(ierr);
514225daaf8SStefano Zampini       ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr);
515225daaf8SStefano Zampini       ierr = PetscContainerSetUserDestroy(c,hypre_array_destroy);CHKERRQ(ierr);
516225daaf8SStefano Zampini       ierr = PetscObjectCompose((PetscObject)(*B),names[i],(PetscObject)c);CHKERRQ(ierr);
517225daaf8SStefano Zampini       ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
518225daaf8SStefano Zampini     }
51963c07aadSStefano Zampini   }
52063c07aadSStefano Zampini   PetscFunctionReturn(0);
52163c07aadSStefano Zampini }
52263c07aadSStefano Zampini 
523613e5ff0Sstefano_zampini static PetscErrorCode MatAIJGetParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
524c1a070e6SStefano Zampini {
525613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *tA;
526c1a070e6SStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
527c1a070e6SStefano Zampini   Mat_SeqAIJ         *diag,*offd;
528c1a070e6SStefano Zampini   PetscInt           *garray,noffd,dnnz,onnz,*row_starts,*col_starts;
529c1a070e6SStefano Zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
530613e5ff0Sstefano_zampini   PetscBool          ismpiaij,isseqaij;
531c1a070e6SStefano Zampini   PetscErrorCode     ierr;
532c1a070e6SStefano Zampini 
533c1a070e6SStefano Zampini   PetscFunctionBegin;
534c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
5354099cc6bSBarry Smith   ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
536c1a070e6SStefano Zampini   if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type);
537c1a070e6SStefano Zampini   if (ismpiaij) {
538c1a070e6SStefano Zampini     Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data);
539c1a070e6SStefano Zampini 
540c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)a->A->data;
541c1a070e6SStefano Zampini     offd   = (Mat_SeqAIJ*)a->B->data;
542c1a070e6SStefano Zampini     garray = a->garray;
543c1a070e6SStefano Zampini     noffd  = a->B->cmap->N;
544c1a070e6SStefano Zampini     dnnz   = diag->nz;
545c1a070e6SStefano Zampini     onnz   = offd->nz;
546c1a070e6SStefano Zampini   } else {
547c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)A->data;
548c1a070e6SStefano Zampini     offd   = NULL;
549c1a070e6SStefano Zampini     garray = NULL;
550c1a070e6SStefano Zampini     noffd  = 0;
551c1a070e6SStefano Zampini     dnnz   = diag->nz;
552c1a070e6SStefano Zampini     onnz   = 0;
553c1a070e6SStefano Zampini   }
554225daaf8SStefano Zampini 
555c1a070e6SStefano Zampini   /* create a temporary ParCSR */
556c1a070e6SStefano Zampini   if (HYPRE_AssumedPartitionCheck()) {
557c1a070e6SStefano Zampini     PetscMPIInt myid;
558c1a070e6SStefano Zampini 
559c1a070e6SStefano Zampini     ierr       = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
560c1a070e6SStefano Zampini     row_starts = A->rmap->range + myid;
561c1a070e6SStefano Zampini     col_starts = A->cmap->range + myid;
562c1a070e6SStefano Zampini   } else {
563c1a070e6SStefano Zampini     row_starts = A->rmap->range;
564c1a070e6SStefano Zampini     col_starts = A->cmap->range;
565c1a070e6SStefano Zampini   }
5667d968826Sstefano_zampini   tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,(HYPRE_Int*)row_starts,(HYPRE_Int*)col_starts,noffd,dnnz,onnz);
567c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(tA,0);
568c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(tA,0);
569c1a070e6SStefano Zampini 
570225daaf8SStefano Zampini   /* set diagonal part */
571c1a070e6SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(tA);
5727d968826Sstefano_zampini   hypre_CSRMatrixI(hdiag)           = (HYPRE_Int*)diag->i;
5737d968826Sstefano_zampini   hypre_CSRMatrixJ(hdiag)           = (HYPRE_Int*)diag->j;
574c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = diag->a;
575c1a070e6SStefano Zampini   hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz;
576c1a070e6SStefano Zampini   hypre_CSRMatrixSetRownnz(hdiag);
577c1a070e6SStefano Zampini   hypre_CSRMatrixSetDataOwner(hdiag,0);
578c1a070e6SStefano Zampini 
579225daaf8SStefano Zampini   /* set offdiagonal part */
580c1a070e6SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(tA);
581c1a070e6SStefano Zampini   if (offd) {
5827d968826Sstefano_zampini     hypre_CSRMatrixI(hoffd)           = (HYPRE_Int*)offd->i;
5837d968826Sstefano_zampini     hypre_CSRMatrixJ(hoffd)           = (HYPRE_Int*)offd->j;
584c1a070e6SStefano Zampini     hypre_CSRMatrixData(hoffd)        = offd->a;
585c1a070e6SStefano Zampini     hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz;
586c1a070e6SStefano Zampini     hypre_CSRMatrixSetRownnz(hoffd);
587c1a070e6SStefano Zampini     hypre_CSRMatrixSetDataOwner(hoffd,0);
588c1a070e6SStefano Zampini     hypre_ParCSRMatrixSetNumNonzeros(tA);
5897d968826Sstefano_zampini     hypre_ParCSRMatrixColMapOffd(tA) = (HYPRE_Int*)garray;
590c1a070e6SStefano Zampini   }
591613e5ff0Sstefano_zampini   *hA = tA;
592613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
593613e5ff0Sstefano_zampini }
594c1a070e6SStefano Zampini 
595613e5ff0Sstefano_zampini static PetscErrorCode MatAIJRestoreParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
596613e5ff0Sstefano_zampini {
597613e5ff0Sstefano_zampini   hypre_CSRMatrix    *hdiag,*hoffd;
598c1a070e6SStefano Zampini 
599613e5ff0Sstefano_zampini   PetscFunctionBegin;
600613e5ff0Sstefano_zampini   hdiag = hypre_ParCSRMatrixDiag(*hA);
601613e5ff0Sstefano_zampini   hoffd = hypre_ParCSRMatrixOffd(*hA);
602c1a070e6SStefano Zampini   /* set pointers to NULL before destroying tA */
603c1a070e6SStefano Zampini   hypre_CSRMatrixI(hdiag)           = NULL;
604c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hdiag)           = NULL;
605c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = NULL;
606c1a070e6SStefano Zampini   hypre_CSRMatrixI(hoffd)           = NULL;
607c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hoffd)           = NULL;
608c1a070e6SStefano Zampini   hypre_CSRMatrixData(hoffd)        = NULL;
609613e5ff0Sstefano_zampini   hypre_ParCSRMatrixColMapOffd(*hA) = NULL;
610613e5ff0Sstefano_zampini   hypre_ParCSRMatrixDestroy(*hA);
611613e5ff0Sstefano_zampini   *hA = NULL;
612613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
613613e5ff0Sstefano_zampini }
614613e5ff0Sstefano_zampini 
615613e5ff0Sstefano_zampini /* calls RAP from BoomerAMG:
6163dad0653Sstefano_zampini    the resulting ParCSR will not own the column and row starts
6173dad0653Sstefano_zampini    It looks like we don't need to have the diagonal entries
6183dad0653Sstefano_zampini    ordered first in the rows of the diagonal part
6193dad0653Sstefano_zampini    for boomerAMGBuildCoarseOperator to work */
620a055b5aaSBarry Smith static PetscErrorCode MatHYPRE_ParCSR_RAP(hypre_ParCSRMatrix *hR, hypre_ParCSRMatrix *hA,hypre_ParCSRMatrix *hP, hypre_ParCSRMatrix **hRAP)
621613e5ff0Sstefano_zampini {
622613e5ff0Sstefano_zampini   PetscErrorCode ierr;
623613e5ff0Sstefano_zampini   HYPRE_Int      P_owns_col_starts,R_owns_row_starts;
624613e5ff0Sstefano_zampini 
625613e5ff0Sstefano_zampini   PetscFunctionBegin;
626613e5ff0Sstefano_zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(hP);
627613e5ff0Sstefano_zampini   R_owns_row_starts = hypre_ParCSRMatrixOwnsRowStarts(hR);
628613e5ff0Sstefano_zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(hR,hA,hP,hRAP));
629613e5ff0Sstefano_zampini   PetscStackCallStandard(hypre_ParCSRMatrixSetNumNonzeros,(*hRAP));
630613e5ff0Sstefano_zampini   /* hypre_BoomerAMGBuildCoarseOperator steals the col_starts from P and the row_starts from R */
631613e5ff0Sstefano_zampini   hypre_ParCSRMatrixSetRowStartsOwner(*hRAP,0);
632613e5ff0Sstefano_zampini   hypre_ParCSRMatrixSetColStartsOwner(*hRAP,0);
633613e5ff0Sstefano_zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(hP,1);
634613e5ff0Sstefano_zampini   if (R_owns_row_starts) hypre_ParCSRMatrixSetRowStartsOwner(hR,1);
635613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
636613e5ff0Sstefano_zampini }
637613e5ff0Sstefano_zampini 
6386f231fbdSstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat P,Mat C)
639613e5ff0Sstefano_zampini {
6406f231fbdSstefano_zampini   Mat                B;
641613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *hA,*hP,*hPtAP;
642613e5ff0Sstefano_zampini   PetscErrorCode     ierr;
643613e5ff0Sstefano_zampini 
644613e5ff0Sstefano_zampini   PetscFunctionBegin;
645613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
646613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(P,&hP);CHKERRQ(ierr);
647613e5ff0Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(hP,hA,hP,&hPtAP);CHKERRQ(ierr);
6486f231fbdSstefano_zampini   ierr = MatCreateFromParCSR(hPtAP,MATAIJ,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
6496f231fbdSstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
650613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
651613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(P,&hP);CHKERRQ(ierr);
6526f231fbdSstefano_zampini   PetscFunctionReturn(0);
6536f231fbdSstefano_zampini }
6546f231fbdSstefano_zampini 
6556f231fbdSstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat P,PetscReal fill,Mat *C)
6566f231fbdSstefano_zampini {
6576f231fbdSstefano_zampini   PetscErrorCode ierr;
658ab4d48faSStefano Zampini 
6596f231fbdSstefano_zampini   PetscFunctionBegin;
6606f231fbdSstefano_zampini   ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
6616f231fbdSstefano_zampini   ierr                   = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
6626f231fbdSstefano_zampini   (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_AIJ_wHYPRE;
663613e5ff0Sstefano_zampini   PetscFunctionReturn(0);
664613e5ff0Sstefano_zampini }
665613e5ff0Sstefano_zampini 
6664cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_AIJ_HYPRE(Mat A,Mat P,Mat C)
667613e5ff0Sstefano_zampini {
6684cc28894Sstefano_zampini   Mat                B;
6694cc28894Sstefano_zampini   Mat_HYPRE          *hP;
670613e5ff0Sstefano_zampini   hypre_ParCSRMatrix *hA = NULL,*Pparcsr,*ptapparcsr;
671613e5ff0Sstefano_zampini   HYPRE_Int          type;
672613e5ff0Sstefano_zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
6734cc28894Sstefano_zampini   PetscBool          ishypre;
674613e5ff0Sstefano_zampini   PetscErrorCode     ierr;
675613e5ff0Sstefano_zampini 
676613e5ff0Sstefano_zampini   PetscFunctionBegin;
6774cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
6784cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(comm,PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
6794cc28894Sstefano_zampini   hP = (Mat_HYPRE*)P->data;
680613e5ff0Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
681613e5ff0Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
682613e5ff0Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
683613e5ff0Sstefano_zampini 
684613e5ff0Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
685613e5ff0Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,hA,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
686613e5ff0Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
687225daaf8SStefano Zampini 
6884cc28894Sstefano_zampini   /* create temporary matrix and merge to C */
6894cc28894Sstefano_zampini   ierr = MatCreateFromParCSR(ptapparcsr,((PetscObject)C)->type_name,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
6904cc28894Sstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
6914cc28894Sstefano_zampini   PetscFunctionReturn(0);
6924cc28894Sstefano_zampini }
6934cc28894Sstefano_zampini 
6944cc28894Sstefano_zampini static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
6954cc28894Sstefano_zampini {
6964cc28894Sstefano_zampini   PetscErrorCode ierr;
6974cc28894Sstefano_zampini 
6984cc28894Sstefano_zampini   PetscFunctionBegin;
6994cc28894Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
7004cc28894Sstefano_zampini     const char *deft = MATAIJ;
7014cc28894Sstefano_zampini     char       type[256];
7024cc28894Sstefano_zampini     PetscBool  flg;
7034cc28894Sstefano_zampini 
7044cc28894Sstefano_zampini     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
7054cc28894Sstefano_zampini     ierr = PetscOptionsFList("-matptap_hypre_outtype","Matrix type",NULL,MatList,deft,type,256,&flg);CHKERRQ(ierr);
7064cc28894Sstefano_zampini     ierr = PetscOptionsEnd();CHKERRQ(ierr);
7074cc28894Sstefano_zampini     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7084cc28894Sstefano_zampini     ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
7094cc28894Sstefano_zampini     if (flg) {
7104cc28894Sstefano_zampini       ierr = MatSetType(*C,type);CHKERRQ(ierr);
7114cc28894Sstefano_zampini     } else {
7124cc28894Sstefano_zampini       ierr = MatSetType(*C,deft);CHKERRQ(ierr);
7134cc28894Sstefano_zampini     }
7144cc28894Sstefano_zampini     (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_HYPRE;
7154cc28894Sstefano_zampini     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7164cc28894Sstefano_zampini   }
7174cc28894Sstefano_zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
7184cc28894Sstefano_zampini   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
719613e5ff0Sstefano_zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
720c1a070e6SStefano Zampini   PetscFunctionReturn(0);
721c1a070e6SStefano Zampini }
722c1a070e6SStefano Zampini 
7234cc28894Sstefano_zampini static PetscErrorCode MatPtAPNumeric_HYPRE_HYPRE(Mat A,Mat P,Mat C)
7244cc28894Sstefano_zampini {
7254cc28894Sstefano_zampini   Mat                B;
7264cc28894Sstefano_zampini   hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr;
7274cc28894Sstefano_zampini   Mat_HYPRE          *hA,*hP;
7284cc28894Sstefano_zampini   PetscBool          ishypre;
7294cc28894Sstefano_zampini   HYPRE_Int          type;
7304cc28894Sstefano_zampini   PetscErrorCode     ierr;
7314cc28894Sstefano_zampini 
7324cc28894Sstefano_zampini   PetscFunctionBegin;
7334cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
7344cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)P),PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
7354cc28894Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr);
7364cc28894Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE);
7374cc28894Sstefano_zampini   hA = (Mat_HYPRE*)A->data;
7384cc28894Sstefano_zampini   hP = (Mat_HYPRE*)P->data;
7394cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
7404cc28894Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
7414cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
7424cc28894Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
7434cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
7444cc28894Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
7454cc28894Sstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
7464cc28894Sstefano_zampini   ierr = MatCreateFromParCSR(ptapparcsr,MATHYPRE,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
7474cc28894Sstefano_zampini   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
7484cc28894Sstefano_zampini   PetscFunctionReturn(0);
7494cc28894Sstefano_zampini }
7504cc28894Sstefano_zampini 
751cd8bc7baSStefano Zampini static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
752cd8bc7baSStefano Zampini {
753cd8bc7baSStefano Zampini   PetscErrorCode ierr;
754cd8bc7baSStefano Zampini 
755cd8bc7baSStefano Zampini   PetscFunctionBegin;
7564cc28894Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
7574cc28894Sstefano_zampini     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7584cc28894Sstefano_zampini     ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
7594cc28894Sstefano_zampini     ierr                   = MatSetType(*C,MATHYPRE);CHKERRQ(ierr);
7604cc28894Sstefano_zampini     (*C)->ops->ptapnumeric = MatPtAPNumeric_HYPRE_HYPRE;
7614cc28894Sstefano_zampini     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
7624cc28894Sstefano_zampini   }
763613e5ff0Sstefano_zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
7644cc28894Sstefano_zampini   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
765613e5ff0Sstefano_zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
766cd8bc7baSStefano Zampini   PetscFunctionReturn(0);
767cd8bc7baSStefano Zampini }
768cd8bc7baSStefano Zampini 
769d501dc42Sstefano_zampini /* calls hypre_ParMatmul
770d501dc42Sstefano_zampini    hypre_ParMatMul uses hypre_ParMatrixCreate with the communicator of hA
7713dad0653Sstefano_zampini    hypre_ParMatrixCreate does not duplicate the communicator
7723dad0653Sstefano_zampini    It looks like we don't need to have the diagonal entries
7733dad0653Sstefano_zampini    ordered first in the rows of the diagonal part
7743dad0653Sstefano_zampini    for boomerAMGBuildCoarseOperator to work */
775d501dc42Sstefano_zampini static PetscErrorCode MatHYPRE_ParCSR_MatMatMult(hypre_ParCSRMatrix *hA, hypre_ParCSRMatrix *hB, hypre_ParCSRMatrix **hAB)
776d501dc42Sstefano_zampini {
777d501dc42Sstefano_zampini   PetscFunctionBegin;
778d501dc42Sstefano_zampini   PetscStackPush("hypre_ParMatmul");
779d501dc42Sstefano_zampini   *hAB = hypre_ParMatmul(hA,hB);
780d501dc42Sstefano_zampini   PetscStackPop;
781d501dc42Sstefano_zampini   PetscFunctionReturn(0);
782d501dc42Sstefano_zampini }
783d501dc42Sstefano_zampini 
7845e5acdf2Sstefano_zampini static PetscErrorCode MatMatMultNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C)
7855e5acdf2Sstefano_zampini {
7865e5acdf2Sstefano_zampini   Mat                D;
787d501dc42Sstefano_zampini   hypre_ParCSRMatrix *hA,*hB,*hAB = NULL;
7885e5acdf2Sstefano_zampini   PetscErrorCode     ierr;
7895e5acdf2Sstefano_zampini 
7905e5acdf2Sstefano_zampini   PetscFunctionBegin;
7915e5acdf2Sstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
7925e5acdf2Sstefano_zampini   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
793d501dc42Sstefano_zampini   ierr = MatHYPRE_ParCSR_MatMatMult(hA,hB,&hAB);CHKERRQ(ierr);
7945e5acdf2Sstefano_zampini   ierr = MatCreateFromParCSR(hAB,MATAIJ,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
7955e5acdf2Sstefano_zampini   ierr = MatHeaderMerge(C,&D);CHKERRQ(ierr);
7965e5acdf2Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
7975e5acdf2Sstefano_zampini   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
7985e5acdf2Sstefano_zampini   PetscFunctionReturn(0);
7995e5acdf2Sstefano_zampini }
8005e5acdf2Sstefano_zampini 
8015e5acdf2Sstefano_zampini PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat B,PetscReal fill,Mat *C)
8025e5acdf2Sstefano_zampini {
8035e5acdf2Sstefano_zampini   PetscErrorCode ierr;
80420e1dc0dSstefano_zampini 
8055e5acdf2Sstefano_zampini   PetscFunctionBegin;
8065e5acdf2Sstefano_zampini   ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
8075e5acdf2Sstefano_zampini   ierr                      = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
8085e5acdf2Sstefano_zampini   (*C)->ops->matmultnumeric = MatMatMultNumeric_AIJ_AIJ_wHYPRE;
8095e5acdf2Sstefano_zampini   PetscFunctionReturn(0);
8105e5acdf2Sstefano_zampini }
8115e5acdf2Sstefano_zampini 
812d501dc42Sstefano_zampini static PetscErrorCode MatMatMultNumeric_HYPRE_HYPRE(Mat A,Mat B,Mat C)
813d501dc42Sstefano_zampini {
814d501dc42Sstefano_zampini   Mat                D;
815d501dc42Sstefano_zampini   hypre_ParCSRMatrix *Aparcsr,*Bparcsr,*ABparcsr = NULL;
816d501dc42Sstefano_zampini   Mat_HYPRE          *hA,*hB;
817d501dc42Sstefano_zampini   PetscBool          ishypre;
818d501dc42Sstefano_zampini   HYPRE_Int          type;
819d501dc42Sstefano_zampini   PetscErrorCode     ierr;
820d501dc42Sstefano_zampini 
821d501dc42Sstefano_zampini   PetscFunctionBegin;
822d501dc42Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)B,MATHYPRE,&ishypre);CHKERRQ(ierr);
823d501dc42Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_USER,"B should be of type %s",MATHYPRE);
824d501dc42Sstefano_zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr);
825d501dc42Sstefano_zampini   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE);
826d501dc42Sstefano_zampini   hA = (Mat_HYPRE*)A->data;
827d501dc42Sstefano_zampini   hB = (Mat_HYPRE*)B->data;
828d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
829d501dc42Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
830d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hB->ij,&type));
831d501dc42Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
832d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
833d501dc42Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hB->ij,(void**)&Bparcsr));
834d501dc42Sstefano_zampini   ierr = MatHYPRE_ParCSR_MatMatMult(Aparcsr,Bparcsr,&ABparcsr);CHKERRQ(ierr);
835d501dc42Sstefano_zampini   ierr = MatCreateFromParCSR(ABparcsr,MATHYPRE,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
836d501dc42Sstefano_zampini   /* need to use HeaderReplace because HeaderMerge messes up with the communicator */
837d501dc42Sstefano_zampini   ierr = MatHeaderReplace(C,&D);CHKERRQ(ierr);
838d501dc42Sstefano_zampini   C->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE;
839d501dc42Sstefano_zampini   PetscFunctionReturn(0);
840d501dc42Sstefano_zampini }
841d501dc42Sstefano_zampini 
842d501dc42Sstefano_zampini static PetscErrorCode MatMatMult_HYPRE_HYPRE(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
843d501dc42Sstefano_zampini {
844d501dc42Sstefano_zampini   PetscErrorCode ierr;
845d501dc42Sstefano_zampini 
846d501dc42Sstefano_zampini   PetscFunctionBegin;
847d501dc42Sstefano_zampini   if (scall == MAT_INITIAL_MATRIX) {
848d501dc42Sstefano_zampini     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
849d501dc42Sstefano_zampini     ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
850d501dc42Sstefano_zampini     ierr                      = MatSetType(*C,MATHYPRE);CHKERRQ(ierr);
851d501dc42Sstefano_zampini     (*C)->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE;
852d501dc42Sstefano_zampini     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
853d501dc42Sstefano_zampini   }
854d501dc42Sstefano_zampini   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
855d501dc42Sstefano_zampini   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
856d501dc42Sstefano_zampini   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
857d501dc42Sstefano_zampini   PetscFunctionReturn(0);
858d501dc42Sstefano_zampini }
859d501dc42Sstefano_zampini 
8603dad0653Sstefano_zampini PETSC_INTERN PetscErrorCode MatTransposeMatMatMultNumeric_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,Mat D)
86120e1dc0dSstefano_zampini {
86220e1dc0dSstefano_zampini   Mat                E;
86320e1dc0dSstefano_zampini   hypre_ParCSRMatrix *hA,*hB,*hC,*hABC;
86420e1dc0dSstefano_zampini   PetscErrorCode     ierr;
86520e1dc0dSstefano_zampini 
86620e1dc0dSstefano_zampini   PetscFunctionBegin;
86720e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
86820e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
86920e1dc0dSstefano_zampini   ierr = MatAIJGetParCSR_Private(C,&hC);CHKERRQ(ierr);
87020e1dc0dSstefano_zampini   ierr = MatHYPRE_ParCSR_RAP(hA,hB,hC,&hABC);CHKERRQ(ierr);
87120e1dc0dSstefano_zampini   ierr = MatCreateFromParCSR(hABC,MATAIJ,PETSC_OWN_POINTER,&E);CHKERRQ(ierr);
87220e1dc0dSstefano_zampini   ierr = MatHeaderMerge(D,&E);CHKERRQ(ierr);
87320e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
87420e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
87520e1dc0dSstefano_zampini   ierr = MatAIJRestoreParCSR_Private(C,&hC);CHKERRQ(ierr);
87620e1dc0dSstefano_zampini   PetscFunctionReturn(0);
87720e1dc0dSstefano_zampini }
87820e1dc0dSstefano_zampini 
8793dad0653Sstefano_zampini PETSC_INTERN PetscErrorCode MatTransposeMatMatMultSymbolic_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,PetscReal fill,Mat *D)
88020e1dc0dSstefano_zampini {
88120e1dc0dSstefano_zampini   PetscErrorCode ierr;
88220e1dc0dSstefano_zampini 
88320e1dc0dSstefano_zampini   PetscFunctionBegin;
88420e1dc0dSstefano_zampini   ierr = MatCreate(PetscObjectComm((PetscObject)A),D);CHKERRQ(ierr);
88520e1dc0dSstefano_zampini   ierr = MatSetType(*D,MATAIJ);CHKERRQ(ierr);
88620e1dc0dSstefano_zampini   PetscFunctionReturn(0);
88720e1dc0dSstefano_zampini }
88820e1dc0dSstefano_zampini 
889ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
89063c07aadSStefano Zampini {
89163c07aadSStefano Zampini   PetscErrorCode ierr;
89263c07aadSStefano Zampini 
89363c07aadSStefano Zampini   PetscFunctionBegin;
89463c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
89563c07aadSStefano Zampini   PetscFunctionReturn(0);
89663c07aadSStefano Zampini }
89763c07aadSStefano Zampini 
898ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
89963c07aadSStefano Zampini {
90063c07aadSStefano Zampini   PetscErrorCode ierr;
90163c07aadSStefano Zampini 
90263c07aadSStefano Zampini   PetscFunctionBegin;
90363c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
90463c07aadSStefano Zampini   PetscFunctionReturn(0);
90563c07aadSStefano Zampini }
90663c07aadSStefano Zampini 
907ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
90863c07aadSStefano Zampini {
90963c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
91063c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
91163c07aadSStefano Zampini   hypre_ParVector    *hx,*hy;
91263c07aadSStefano Zampini   PetscScalar        *ax,*ay,*sax,*say;
91363c07aadSStefano Zampini   PetscErrorCode     ierr;
91463c07aadSStefano Zampini 
91563c07aadSStefano Zampini   PetscFunctionBegin;
91663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
91763c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
91863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
91963c07aadSStefano Zampini   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
92063c07aadSStefano Zampini   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
92163c07aadSStefano Zampini   if (trans) {
92258968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ay,say);
92358968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax);
92463c07aadSStefano Zampini     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
92558968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,say,ay);
92658968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax);
92763c07aadSStefano Zampini   } else {
92858968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax);
92958968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ay,say);
93063c07aadSStefano Zampini     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
93158968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax);
93258968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,say,ay);
93363c07aadSStefano Zampini   }
93463c07aadSStefano Zampini   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
93563c07aadSStefano Zampini   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
93663c07aadSStefano Zampini   PetscFunctionReturn(0);
93763c07aadSStefano Zampini }
93863c07aadSStefano Zampini 
939ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A)
94063c07aadSStefano Zampini {
94163c07aadSStefano Zampini   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
94263c07aadSStefano Zampini   PetscErrorCode ierr;
94363c07aadSStefano Zampini 
94463c07aadSStefano Zampini   PetscFunctionBegin;
94563c07aadSStefano Zampini   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
94663c07aadSStefano Zampini   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
947978814f1SStefano Zampini   if (hA->ij) {
948978814f1SStefano Zampini     if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL;
949978814f1SStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
950978814f1SStefano Zampini   }
95163c07aadSStefano Zampini   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); }
95263c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
9532df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr);
954c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr);
955c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr);
956d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPRESetPreallocation_C",NULL);CHKERRQ(ierr);
957dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPREGetParCSR_C",NULL);CHKERRQ(ierr);
958*75d48cdbSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_is_hypre_C",NULL);CHKERRQ(ierr);
95963c07aadSStefano Zampini   ierr = PetscFree(A->data);CHKERRQ(ierr);
96063c07aadSStefano Zampini   PetscFunctionReturn(0);
96163c07aadSStefano Zampini }
96263c07aadSStefano Zampini 
963ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A)
96463c07aadSStefano Zampini {
9654ec6421dSstefano_zampini   PetscErrorCode ierr;
9664ec6421dSstefano_zampini 
9674ec6421dSstefano_zampini   PetscFunctionBegin;
9684ec6421dSstefano_zampini   ierr = MatHYPRESetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr);
9694ec6421dSstefano_zampini   PetscFunctionReturn(0);
9704ec6421dSstefano_zampini }
9714ec6421dSstefano_zampini 
9724ec6421dSstefano_zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
9734ec6421dSstefano_zampini {
97463c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
97563c07aadSStefano Zampini   Vec                x,b;
97663c07aadSStefano Zampini   PetscErrorCode     ierr;
97763c07aadSStefano Zampini 
97863c07aadSStefano Zampini   PetscFunctionBegin;
9794ec6421dSstefano_zampini   if (mode == MAT_FLUSH_ASSEMBLY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_FLUSH_ASSEMBLY currently not supported with MATHYPRE");
9804ec6421dSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
98163c07aadSStefano Zampini   if (hA->x) PetscFunctionReturn(0);
98263c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
98363c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
98463c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
98563c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
98663c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
98763c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
98863c07aadSStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
98963c07aadSStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
99063c07aadSStefano Zampini   PetscFunctionReturn(0);
99163c07aadSStefano Zampini }
99263c07aadSStefano Zampini 
993d975228cSstefano_zampini #define MATHYPRE_SCRATCH 2048
994d975228cSstefano_zampini 
995d975228cSstefano_zampini PetscErrorCode MatSetValues_HYPRE(Mat A, PetscInt nr, const PetscInt rows[], PetscInt nc, const PetscInt cols[], const PetscScalar v[], InsertMode ins)
996d975228cSstefano_zampini {
997d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
998d975228cSstefano_zampini   PetscScalar        *vals = (PetscScalar *)v;
999d975228cSstefano_zampini   PetscScalar        sscr[MATHYPRE_SCRATCH];
10007d968826Sstefano_zampini   HYPRE_Int          cscr[2][MATHYPRE_SCRATCH];
10017d968826Sstefano_zampini   HYPRE_Int          i,nzc;
1002d975228cSstefano_zampini   PetscErrorCode     ierr;
1003d975228cSstefano_zampini 
1004d975228cSstefano_zampini   PetscFunctionBegin;
1005d975228cSstefano_zampini   for (i=0,nzc=0;i<nc;i++) {
1006d975228cSstefano_zampini     if (cols[i] >= 0) {
1007d975228cSstefano_zampini       cscr[0][nzc  ] = cols[i];
1008d975228cSstefano_zampini       cscr[1][nzc++] = i;
1009d975228cSstefano_zampini     }
1010d975228cSstefano_zampini   }
1011d975228cSstefano_zampini   if (!nzc) PetscFunctionReturn(0);
1012d975228cSstefano_zampini 
1013d975228cSstefano_zampini   if (ins == ADD_VALUES) {
1014d975228cSstefano_zampini     for (i=0;i<nr;i++) {
1015e3977e59Sstefano_zampini       if (rows[i] >= 0 && nzc) {
1016d975228cSstefano_zampini         PetscInt j;
1017d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
10187d968826Sstefano_zampini         PetscStackCallStandard(HYPRE_IJMatrixAddToValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
1019d975228cSstefano_zampini       }
1020d975228cSstefano_zampini       vals += nc;
1021d975228cSstefano_zampini     }
1022d975228cSstefano_zampini   } else { /* INSERT_VALUES */
1023d975228cSstefano_zampini #if defined(PETSC_USE_DEBUG)
1024d975228cSstefano_zampini     /* Insert values cannot be used to insert offproc entries */
1025d975228cSstefano_zampini     PetscInt rst,ren;
1026d975228cSstefano_zampini     ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
1027d975228cSstefano_zampini     for (i=0;i<nr;i++)
1028d975228cSstefano_zampini       if (rows[i] >= 0 && (rows[i] < rst || rows[i] >= ren)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use INSERT_VALUES for off-proc entries with MatHYPRE. Use ADD_VALUES instead");
1029d975228cSstefano_zampini #endif
1030d975228cSstefano_zampini     for (i=0;i<nr;i++) {
1031e3977e59Sstefano_zampini       if (rows[i] >= 0 && nzc) {
1032d975228cSstefano_zampini         PetscInt j;
1033d975228cSstefano_zampini         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
10347d968826Sstefano_zampini         PetscStackCallStandard(HYPRE_IJMatrixSetValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
1035d975228cSstefano_zampini       }
1036d975228cSstefano_zampini       vals += nc;
1037d975228cSstefano_zampini     }
1038d975228cSstefano_zampini   }
1039d975228cSstefano_zampini   PetscFunctionReturn(0);
1040d975228cSstefano_zampini }
1041d975228cSstefano_zampini 
1042d975228cSstefano_zampini static PetscErrorCode MatHYPRESetPreallocation_HYPRE(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1043d975228cSstefano_zampini {
1044d975228cSstefano_zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
10457d968826Sstefano_zampini   HYPRE_Int          *hdnnz,*honnz;
104606a29025Sstefano_zampini   PetscInt           i,rs,re,cs,ce,bs;
1047d975228cSstefano_zampini   PetscMPIInt        size;
1048d975228cSstefano_zampini   PetscErrorCode     ierr;
1049d975228cSstefano_zampini 
1050d975228cSstefano_zampini   PetscFunctionBegin;
105106a29025Sstefano_zampini   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1052d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1053d975228cSstefano_zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1054d975228cSstefano_zampini   rs   = A->rmap->rstart;
1055d975228cSstefano_zampini   re   = A->rmap->rend;
1056d975228cSstefano_zampini   cs   = A->cmap->rstart;
1057d975228cSstefano_zampini   ce   = A->cmap->rend;
1058d975228cSstefano_zampini   if (!hA->ij) {
1059d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rs,re-1,cs,ce-1,&hA->ij));
1060d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1061d975228cSstefano_zampini   } else {
10627d968826Sstefano_zampini     HYPRE_Int hrs,hre,hcs,hce;
1063d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixGetLocalRange,(hA->ij,&hrs,&hre,&hcs,&hce));
1064d975228cSstefano_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);
1065d975228cSstefano_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);
1066d975228cSstefano_zampini   }
1067d975228cSstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1068d975228cSstefano_zampini 
106906a29025Sstefano_zampini   if (dnz == PETSC_DEFAULT || dnz == PETSC_DECIDE) dnz = 10*bs;
107006a29025Sstefano_zampini   if (onz == PETSC_DEFAULT || onz == PETSC_DECIDE) onz = 10*bs;
107106a29025Sstefano_zampini 
1072d975228cSstefano_zampini   if (!dnnz) {
1073d975228cSstefano_zampini     ierr = PetscMalloc1(A->rmap->n,&hdnnz);CHKERRQ(ierr);
1074d975228cSstefano_zampini     for (i=0;i<A->rmap->n;i++) hdnnz[i] = dnz;
1075d975228cSstefano_zampini   } else {
10767d968826Sstefano_zampini     hdnnz = (HYPRE_Int*)dnnz;
1077d975228cSstefano_zampini   }
1078d975228cSstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1079d975228cSstefano_zampini   if (size > 1) {
1080d975228cSstefano_zampini     if (!onnz) {
1081d975228cSstefano_zampini       ierr = PetscMalloc1(A->rmap->n,&honnz);CHKERRQ(ierr);
1082d975228cSstefano_zampini       for (i=0;i<A->rmap->n;i++) honnz[i] = onz;
1083d975228cSstefano_zampini     } else {
10847d968826Sstefano_zampini       honnz = (HYPRE_Int*)onnz;
1085d975228cSstefano_zampini     }
1086d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(hA->ij,hdnnz,honnz));
1087d975228cSstefano_zampini   } else {
1088d975228cSstefano_zampini     honnz = NULL;
1089d975228cSstefano_zampini     PetscStackCallStandard(HYPRE_IJMatrixSetRowSizes,(hA->ij,hdnnz));
1090d975228cSstefano_zampini   }
1091d975228cSstefano_zampini   if (!dnnz) {
1092d975228cSstefano_zampini     ierr = PetscFree(hdnnz);CHKERRQ(ierr);
1093d975228cSstefano_zampini   }
1094d975228cSstefano_zampini   if (!onnz && honnz) {
1095d975228cSstefano_zampini     ierr = PetscFree(honnz);CHKERRQ(ierr);
1096d975228cSstefano_zampini   }
109706a29025Sstefano_zampini   A->preallocated = PETSC_TRUE;
1098d975228cSstefano_zampini 
1099d975228cSstefano_zampini   /* SetDiagOffdSizes sets hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0 */
1100d975228cSstefano_zampini   {
1101d975228cSstefano_zampini     hypre_AuxParCSRMatrix *aux_matrix;
1102d975228cSstefano_zampini     aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij);
1103d975228cSstefano_zampini     hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1;
1104d975228cSstefano_zampini   }
1105d975228cSstefano_zampini   PetscFunctionReturn(0);
1106d975228cSstefano_zampini }
1107d975228cSstefano_zampini 
1108d975228cSstefano_zampini /*@C
1109d975228cSstefano_zampini    MatHYPRESetPreallocation - Preallocates memory for a sparse parallel matrix in HYPRE IJ format
1110d975228cSstefano_zampini 
1111d975228cSstefano_zampini    Collective on Mat
1112d975228cSstefano_zampini 
1113d975228cSstefano_zampini    Input Parameters:
1114d975228cSstefano_zampini +  A - the matrix
1115d975228cSstefano_zampini .  dnz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1116d975228cSstefano_zampini           (same value is used for all local rows)
1117d975228cSstefano_zampini .  dnnz - array containing the number of nonzeros in the various rows of the
1118d975228cSstefano_zampini           DIAGONAL portion of the local submatrix (possibly different for each row)
1119d975228cSstefano_zampini           or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure.
1120d975228cSstefano_zampini           The size of this array is equal to the number of local rows, i.e 'm'.
1121d975228cSstefano_zampini           For matrices that will be factored, you must leave room for (and set)
1122d975228cSstefano_zampini           the diagonal entry even if it is zero.
1123d975228cSstefano_zampini .  onz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1124d975228cSstefano_zampini           submatrix (same value is used for all local rows).
1125d975228cSstefano_zampini -  onnz - array containing the number of nonzeros in the various rows of the
1126d975228cSstefano_zampini           OFF-DIAGONAL portion of the local submatrix (possibly different for
1127d975228cSstefano_zampini           each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero
1128d975228cSstefano_zampini           structure. The size of this array is equal to the number
1129d975228cSstefano_zampini           of local rows, i.e 'm'.
1130d975228cSstefano_zampini 
113195452b02SPatrick Sanan    Notes:
113295452b02SPatrick Sanan     If the *nnz parameter is given then the *nz parameter is ignored; for sequential matrices, onz and onnz are ignored.
1133d975228cSstefano_zampini 
1134d975228cSstefano_zampini    Level: intermediate
1135d975228cSstefano_zampini 
1136d975228cSstefano_zampini .keywords: matrix, aij, compressed row, sparse, parallel
1137d975228cSstefano_zampini 
1138d975228cSstefano_zampini .seealso: MatCreate(), MatMPIAIJSetPreallocation, MATHYPRE
1139d975228cSstefano_zampini @*/
1140d975228cSstefano_zampini PetscErrorCode MatHYPRESetPreallocation(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1141d975228cSstefano_zampini {
1142d975228cSstefano_zampini   PetscErrorCode ierr;
1143d975228cSstefano_zampini 
1144d975228cSstefano_zampini   PetscFunctionBegin;
1145d975228cSstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1146d975228cSstefano_zampini   PetscValidType(A,1);
1147d975228cSstefano_zampini   ierr = PetscTryMethod(A,"MatHYPRESetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(A,dnz,dnnz,onz,onnz));CHKERRQ(ierr);
1148d975228cSstefano_zampini   PetscFunctionReturn(0);
1149d975228cSstefano_zampini }
1150d975228cSstefano_zampini 
1151225daaf8SStefano Zampini /*
1152225daaf8SStefano Zampini    MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix
1153225daaf8SStefano Zampini 
1154225daaf8SStefano Zampini    Collective
1155225daaf8SStefano Zampini 
1156225daaf8SStefano Zampini    Input Parameters:
1157225daaf8SStefano Zampini +  vparcsr  - the pointer to the hypre_ParCSRMatrix
1158bb4689ddSStefano Zampini .  mtype    - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported.
1159225daaf8SStefano Zampini -  copymode - PETSc copying options
1160225daaf8SStefano Zampini 
1161225daaf8SStefano Zampini    Output Parameter:
1162225daaf8SStefano Zampini .  A  - the matrix
1163225daaf8SStefano Zampini 
1164225daaf8SStefano Zampini    Level: intermediate
1165225daaf8SStefano Zampini 
1166225daaf8SStefano Zampini .seealso: MatHYPRE, PetscCopyMode
1167225daaf8SStefano Zampini */
1168dd9c0a25Sstefano_zampini PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(hypre_ParCSRMatrix *vparcsr, MatType mtype, PetscCopyMode copymode, Mat* A)
1169978814f1SStefano Zampini {
1170225daaf8SStefano Zampini   Mat                   T;
1171978814f1SStefano Zampini   Mat_HYPRE             *hA;
117263c07aadSStefano Zampini   hypre_ParCSRMatrix    *parcsr;
1173978814f1SStefano Zampini   MPI_Comm              comm;
1174978814f1SStefano Zampini   PetscInt              rstart,rend,cstart,cend,M,N;
1175bb4689ddSStefano Zampini   PetscBool             isseqaij,ismpiaij,isaij,ishyp,isis;
1176978814f1SStefano Zampini   PetscErrorCode        ierr;
1177978814f1SStefano Zampini 
1178978814f1SStefano Zampini   PetscFunctionBegin;
1179978814f1SStefano Zampini   parcsr = (hypre_ParCSRMatrix *)vparcsr;
1180978814f1SStefano Zampini   comm   = hypre_ParCSRMatrixComm(parcsr);
1181225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
1182225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
1183225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr);
1184225daaf8SStefano Zampini   ierr   = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr);
11858cfe8d00SStefano Zampini   ierr   = PetscStrcmp(mtype,MATIS,&isis);CHKERRQ(ierr);
1186225daaf8SStefano Zampini   isaij  = (PetscBool)(isseqaij || ismpiaij || isaij);
1187bb4689ddSStefano 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);
1188225daaf8SStefano Zampini   if (ishyp && copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES");
1189978814f1SStefano Zampini 
1190978814f1SStefano Zampini   /* access ParCSRMatrix */
1191978814f1SStefano Zampini   rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr);
1192978814f1SStefano Zampini   rend   = hypre_ParCSRMatrixLastRowIndex(parcsr);
1193978814f1SStefano Zampini   cstart = hypre_ParCSRMatrixFirstColDiag(parcsr);
1194978814f1SStefano Zampini   cend   = hypre_ParCSRMatrixLastColDiag(parcsr);
1195978814f1SStefano Zampini   M      = hypre_ParCSRMatrixGlobalNumRows(parcsr);
1196978814f1SStefano Zampini   N      = hypre_ParCSRMatrixGlobalNumCols(parcsr);
1197978814f1SStefano Zampini 
1198fa92c42cSstefano_zampini   /* fix for empty local rows/columns */
1199fa92c42cSstefano_zampini   if (rend < rstart) rend = rstart;
1200fa92c42cSstefano_zampini   if (cend < cstart) cend = cstart;
1201fa92c42cSstefano_zampini 
1202978814f1SStefano Zampini   /* create PETSc matrix with MatHYPRE */
1203225daaf8SStefano Zampini   ierr = MatCreate(comm,&T);CHKERRQ(ierr);
1204225daaf8SStefano Zampini   ierr = MatSetSizes(T,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr);
1205225daaf8SStefano Zampini   ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr);
1206225daaf8SStefano Zampini   hA   = (Mat_HYPRE*)(T->data);
1207978814f1SStefano Zampini 
1208978814f1SStefano Zampini   /* create HYPRE_IJMatrix */
1209978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
1210978814f1SStefano Zampini 
1211978814f1SStefano Zampini   /* set ParCSR object */
1212978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1213978814f1SStefano Zampini   hypre_IJMatrixObject(hA->ij) = parcsr;
12144ec6421dSstefano_zampini   T->preallocated = PETSC_TRUE;
1215978814f1SStefano Zampini 
1216978814f1SStefano Zampini   /* set assembled flag */
1217978814f1SStefano Zampini   hypre_IJMatrixAssembleFlag(hA->ij) = 1;
1218978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1219225daaf8SStefano Zampini   if (ishyp) {
12206d2a658fSstefano_zampini     PetscMPIInt myid = 0;
12216d2a658fSstefano_zampini 
12226d2a658fSstefano_zampini     /* make sure we always have row_starts and col_starts available */
12236d2a658fSstefano_zampini     if (HYPRE_AssumedPartitionCheck()) {
12246d2a658fSstefano_zampini       ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
12256d2a658fSstefano_zampini     }
12266d2a658fSstefano_zampini     if (!hypre_ParCSRMatrixOwnsColStarts(parcsr)) {
12276d2a658fSstefano_zampini       PetscLayout map;
12286d2a658fSstefano_zampini 
12296d2a658fSstefano_zampini       ierr = MatGetLayouts(T,NULL,&map);CHKERRQ(ierr);
12306d2a658fSstefano_zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
12311c92d2d0Sstefano_zampini       hypre_ParCSRMatrixColStarts(parcsr) = (HYPRE_Int*)(map->range + myid);
12326d2a658fSstefano_zampini     }
12336d2a658fSstefano_zampini     if (!hypre_ParCSRMatrixOwnsRowStarts(parcsr)) {
12346d2a658fSstefano_zampini       PetscLayout map;
12356d2a658fSstefano_zampini 
12366d2a658fSstefano_zampini       ierr = MatGetLayouts(T,&map,NULL);CHKERRQ(ierr);
12376d2a658fSstefano_zampini       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
12381c92d2d0Sstefano_zampini       hypre_ParCSRMatrixRowStarts(parcsr) = (HYPRE_Int*)(map->range + myid);
12396d2a658fSstefano_zampini     }
1240978814f1SStefano Zampini     /* prevent from freeing the pointer */
1241978814f1SStefano Zampini     if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE;
1242225daaf8SStefano Zampini     *A   = T;
1243978814f1SStefano Zampini     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1244978814f1SStefano Zampini     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1245bb4689ddSStefano Zampini   } else if (isaij) {
1246bb4689ddSStefano Zampini     if (copymode != PETSC_OWN_POINTER) {
1247225daaf8SStefano Zampini       /* prevent from freeing the pointer */
1248225daaf8SStefano Zampini       hA->inner_free = PETSC_FALSE;
1249225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
1250225daaf8SStefano Zampini       ierr = MatDestroy(&T);CHKERRQ(ierr);
1251225daaf8SStefano Zampini     } else { /* AIJ return type with PETSC_OWN_POINTER */
1252225daaf8SStefano Zampini       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1253225daaf8SStefano Zampini       *A   = T;
1254225daaf8SStefano Zampini     }
1255bb4689ddSStefano Zampini   } else if (isis) {
12568cfe8d00SStefano Zampini     ierr = MatConvert_HYPRE_IS(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
12578cfe8d00SStefano Zampini     if (copymode != PETSC_OWN_POINTER) hA->inner_free = PETSC_FALSE;
1258bb4689ddSStefano Zampini     ierr = MatDestroy(&T);CHKERRQ(ierr);
1259bb4689ddSStefano Zampini   }
1260978814f1SStefano Zampini   PetscFunctionReturn(0);
1261978814f1SStefano Zampini }
1262978814f1SStefano Zampini 
126368ec7858SStefano Zampini static PetscErrorCode MatHYPREGetParCSR_HYPRE(Mat A, hypre_ParCSRMatrix **parcsr)
1264dd9c0a25Sstefano_zampini {
1265dd9c0a25Sstefano_zampini   Mat_HYPRE*            hA = (Mat_HYPRE*)A->data;
1266dd9c0a25Sstefano_zampini   HYPRE_Int             type;
1267dd9c0a25Sstefano_zampini   PetscErrorCode        ierr;
1268dd9c0a25Sstefano_zampini 
1269dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1270dd9c0a25Sstefano_zampini   if (!hA->ij) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"HYPRE_IJMatrix not present");
1271dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
1272dd9c0a25Sstefano_zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"HYPRE_IJMatrix is not of type HYPRE_PARCSR");
1273dd9c0a25Sstefano_zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)parcsr));
1274dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1275dd9c0a25Sstefano_zampini }
1276dd9c0a25Sstefano_zampini 
1277dd9c0a25Sstefano_zampini /*
1278dd9c0a25Sstefano_zampini    MatHYPREGetParCSR - Gets the pointer to the ParCSR matrix
1279dd9c0a25Sstefano_zampini 
1280dd9c0a25Sstefano_zampini    Not collective
1281dd9c0a25Sstefano_zampini 
1282dd9c0a25Sstefano_zampini    Input Parameters:
1283dd9c0a25Sstefano_zampini +  A  - the MATHYPRE object
1284dd9c0a25Sstefano_zampini 
1285dd9c0a25Sstefano_zampini    Output Parameter:
1286dd9c0a25Sstefano_zampini .  parcsr  - the pointer to the hypre_ParCSRMatrix
1287dd9c0a25Sstefano_zampini 
1288dd9c0a25Sstefano_zampini    Level: intermediate
1289dd9c0a25Sstefano_zampini 
1290dd9c0a25Sstefano_zampini .seealso: MatHYPRE, PetscCopyMode
1291dd9c0a25Sstefano_zampini */
1292dd9c0a25Sstefano_zampini PetscErrorCode MatHYPREGetParCSR(Mat A, hypre_ParCSRMatrix **parcsr)
1293dd9c0a25Sstefano_zampini {
1294dd9c0a25Sstefano_zampini   PetscErrorCode ierr;
1295dd9c0a25Sstefano_zampini 
1296dd9c0a25Sstefano_zampini   PetscFunctionBegin;
1297dd9c0a25Sstefano_zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1298dd9c0a25Sstefano_zampini   PetscValidType(A,1);
1299dd9c0a25Sstefano_zampini   ierr = PetscUseMethod(A,"MatHYPREGetParCSR_C",(Mat,hypre_ParCSRMatrix**),(A,parcsr));CHKERRQ(ierr);
1300dd9c0a25Sstefano_zampini   PetscFunctionReturn(0);
1301dd9c0a25Sstefano_zampini }
1302dd9c0a25Sstefano_zampini 
130368ec7858SStefano Zampini static PetscErrorCode MatMissingDiagonal_HYPRE(Mat A, PetscBool *missing, PetscInt *dd)
130468ec7858SStefano Zampini {
130568ec7858SStefano Zampini   hypre_ParCSRMatrix *parcsr;
130668ec7858SStefano Zampini   hypre_CSRMatrix    *ha;
130768ec7858SStefano Zampini   PetscInt           rst;
130868ec7858SStefano Zampini   PetscErrorCode     ierr;
130968ec7858SStefano Zampini 
131068ec7858SStefano Zampini   PetscFunctionBegin;
131168ec7858SStefano Zampini   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented with non-square diagonal blocks");
131268ec7858SStefano Zampini   ierr = MatGetOwnershipRange(A,&rst,NULL);CHKERRQ(ierr);
131368ec7858SStefano Zampini   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
131468ec7858SStefano Zampini   if (missing) *missing = PETSC_FALSE;
131568ec7858SStefano Zampini   if (dd) *dd = -1;
131668ec7858SStefano Zampini   ha = hypre_ParCSRMatrixDiag(parcsr);
131768ec7858SStefano Zampini   if (ha) {
131868299464SStefano Zampini     PetscInt  size,i;
131968299464SStefano Zampini     HYPRE_Int *ii,*jj;
132068ec7858SStefano Zampini 
132168ec7858SStefano Zampini     size = hypre_CSRMatrixNumRows(ha);
132268ec7858SStefano Zampini     ii   = hypre_CSRMatrixI(ha);
132368ec7858SStefano Zampini     jj   = hypre_CSRMatrixJ(ha);
132468ec7858SStefano Zampini     for (i = 0; i < size; i++) {
132568ec7858SStefano Zampini       PetscInt  j;
132668ec7858SStefano Zampini       PetscBool found = PETSC_FALSE;
132768ec7858SStefano Zampini 
132868ec7858SStefano Zampini       for (j = ii[i]; j < ii[i+1] && !found; j++)
132968ec7858SStefano Zampini         found = (jj[j] == i) ? PETSC_TRUE : PETSC_FALSE;
133068ec7858SStefano Zampini 
133168ec7858SStefano Zampini       if (!found) {
133268ec7858SStefano Zampini         PetscInfo1(A,"Matrix is missing local diagonal entry %D\n",i);
133368ec7858SStefano Zampini         if (missing) *missing = PETSC_TRUE;
133468ec7858SStefano Zampini         if (dd) *dd = i+rst;
133568ec7858SStefano Zampini         PetscFunctionReturn(0);
133668ec7858SStefano Zampini       }
133768ec7858SStefano Zampini     }
133868ec7858SStefano Zampini     if (!size) {
133968ec7858SStefano Zampini       PetscInfo(A,"Matrix has no diagonal entries therefore is missing diagonal\n");
134068ec7858SStefano Zampini       if (missing) *missing = PETSC_TRUE;
134168ec7858SStefano Zampini       if (dd) *dd = rst;
134268ec7858SStefano Zampini     }
134368ec7858SStefano Zampini   } else {
134468ec7858SStefano Zampini     PetscInfo(A,"Matrix has no diagonal entries therefore is missing diagonal\n");
134568ec7858SStefano Zampini     if (missing) *missing = PETSC_TRUE;
134668ec7858SStefano Zampini     if (dd) *dd = rst;
134768ec7858SStefano Zampini   }
134868ec7858SStefano Zampini   PetscFunctionReturn(0);
134968ec7858SStefano Zampini }
135068ec7858SStefano Zampini 
135168ec7858SStefano Zampini static PetscErrorCode MatScale_HYPRE(Mat A, PetscScalar s)
135268ec7858SStefano Zampini {
135368ec7858SStefano Zampini   hypre_ParCSRMatrix *parcsr;
135468ec7858SStefano Zampini   hypre_CSRMatrix    *ha;
135568ec7858SStefano Zampini   PetscErrorCode     ierr;
135668ec7858SStefano Zampini 
135768ec7858SStefano Zampini   PetscFunctionBegin;
135868ec7858SStefano Zampini   ierr  = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
135968ec7858SStefano Zampini   /* diagonal part */
136068ec7858SStefano Zampini   ha = hypre_ParCSRMatrixDiag(parcsr);
136168ec7858SStefano Zampini   if (ha) {
136268299464SStefano Zampini     PetscInt    size,i;
136368299464SStefano Zampini     HYPRE_Int   *ii;
136468ec7858SStefano Zampini     PetscScalar *a;
136568ec7858SStefano Zampini 
136668ec7858SStefano Zampini     size = hypre_CSRMatrixNumRows(ha);
136768ec7858SStefano Zampini     a    = hypre_CSRMatrixData(ha);
136868ec7858SStefano Zampini     ii   = hypre_CSRMatrixI(ha);
136968ec7858SStefano Zampini     for (i = 0; i < ii[size]; i++) a[i] *= s;
137068ec7858SStefano Zampini   }
137168ec7858SStefano Zampini   /* offdiagonal part */
137268ec7858SStefano Zampini   ha = hypre_ParCSRMatrixOffd(parcsr);
137368ec7858SStefano Zampini   if (ha) {
137468299464SStefano Zampini     PetscInt    size,i;
137568299464SStefano Zampini     HYPRE_Int   *ii;
137668ec7858SStefano Zampini     PetscScalar *a;
137768ec7858SStefano Zampini 
137868ec7858SStefano Zampini     size = hypre_CSRMatrixNumRows(ha);
137968ec7858SStefano Zampini     a    = hypre_CSRMatrixData(ha);
138068ec7858SStefano Zampini     ii   = hypre_CSRMatrixI(ha);
138168ec7858SStefano Zampini     for (i = 0; i < ii[size]; i++) a[i] *= s;
138268ec7858SStefano Zampini   }
138368ec7858SStefano Zampini   PetscFunctionReturn(0);
138468ec7858SStefano Zampini }
138568ec7858SStefano Zampini 
138668ec7858SStefano Zampini static PetscErrorCode MatZeroRowsColumns_HYPRE(Mat A, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
138768ec7858SStefano Zampini {
138868ec7858SStefano Zampini   hypre_ParCSRMatrix *parcsr;
138968299464SStefano Zampini   HYPRE_Int          *lrows;
139068299464SStefano Zampini   PetscInt           rst,ren,i;
139168ec7858SStefano Zampini   PetscErrorCode     ierr;
139268ec7858SStefano Zampini 
139368ec7858SStefano Zampini   PetscFunctionBegin;
139468ec7858SStefano Zampini   if (x || b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");
139568ec7858SStefano Zampini   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
139668ec7858SStefano Zampini   ierr = PetscMalloc1(numRows,&lrows);CHKERRQ(ierr);
139768ec7858SStefano Zampini   ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
139868ec7858SStefano Zampini   for (i=0;i<numRows;i++) {
139968ec7858SStefano Zampini     if (rows[i] < rst || rows[i] >= ren)
140068ec7858SStefano Zampini       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Non-local rows not yet supported");
140168ec7858SStefano Zampini     lrows[i] = rows[i] - rst;
140268ec7858SStefano Zampini   }
140368ec7858SStefano Zampini   PetscStackCallStandard(hypre_ParCSRMatrixEliminateRowsCols,(parcsr,numRows,lrows));
140468ec7858SStefano Zampini   ierr = PetscFree(lrows);CHKERRQ(ierr);
140568ec7858SStefano Zampini   PetscFunctionReturn(0);
140668ec7858SStefano Zampini }
140768ec7858SStefano Zampini 
1408a055b5aaSBarry Smith /*MC
1409a055b5aaSBarry Smith    MATHYPRE - MATHYPRE = "hypre" - A matrix type to be used for sequential and parallel sparse matrices
1410a055b5aaSBarry Smith           based on the hypre IJ interface.
1411a055b5aaSBarry Smith 
1412a055b5aaSBarry Smith    Level: intermediate
1413a055b5aaSBarry Smith 
1414a055b5aaSBarry Smith .seealso: MatCreate()
1415a055b5aaSBarry Smith M*/
1416a055b5aaSBarry Smith 
141763c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
141863c07aadSStefano Zampini {
141963c07aadSStefano Zampini   Mat_HYPRE      *hB;
142063c07aadSStefano Zampini   PetscErrorCode ierr;
142163c07aadSStefano Zampini 
142263c07aadSStefano Zampini   PetscFunctionBegin;
142363c07aadSStefano Zampini   ierr           = PetscNewLog(B,&hB);CHKERRQ(ierr);
1424978814f1SStefano Zampini   hB->inner_free = PETSC_TRUE;
1425978814f1SStefano Zampini 
142663c07aadSStefano Zampini   B->data       = (void*)hB;
142763c07aadSStefano Zampini   B->rmap->bs   = 1;
142863c07aadSStefano Zampini   B->assembled  = PETSC_FALSE;
142963c07aadSStefano Zampini 
1430de393100SStefano Zampini   ierr = PetscMemzero(B->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
143163c07aadSStefano Zampini   B->ops->mult            = MatMult_HYPRE;
143263c07aadSStefano Zampini   B->ops->multtranspose   = MatMultTranspose_HYPRE;
143363c07aadSStefano Zampini   B->ops->setup           = MatSetUp_HYPRE;
143463c07aadSStefano Zampini   B->ops->destroy         = MatDestroy_HYPRE;
143563c07aadSStefano Zampini   B->ops->assemblyend     = MatAssemblyEnd_HYPRE;
1436cd8bc7baSStefano Zampini   B->ops->ptap            = MatPtAP_HYPRE_HYPRE;
1437d501dc42Sstefano_zampini   B->ops->matmult         = MatMatMult_HYPRE_HYPRE;
1438d975228cSstefano_zampini   B->ops->setvalues       = MatSetValues_HYPRE;
143968ec7858SStefano Zampini   B->ops->missingdiagonal = MatMissingDiagonal_HYPRE;
144068ec7858SStefano Zampini   B->ops->scale           = MatScale_HYPRE;
144168ec7858SStefano Zampini   B->ops->zerorowscolumns = MatZeroRowsColumns_HYPRE;
144263c07aadSStefano Zampini 
144363c07aadSStefano Zampini   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
144463c07aadSStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
144563c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
14462df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr);
1447c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1448c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1449d975228cSstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPRESetPreallocation_C",MatHYPRESetPreallocation_HYPRE);CHKERRQ(ierr);
1450dd9c0a25Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPREGetParCSR_C",MatHYPREGetParCSR_HYPRE);CHKERRQ(ierr);
1451*75d48cdbSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_is_hypre_C",MatPtAP_IS_XAIJ);CHKERRQ(ierr);
145263c07aadSStefano Zampini   PetscFunctionReturn(0);
145363c07aadSStefano Zampini }
145463c07aadSStefano Zampini 
1455225daaf8SStefano Zampini static PetscErrorCode hypre_array_destroy(void *ptr)
1456225daaf8SStefano Zampini {
1457225daaf8SStefano Zampini    PetscFunctionBegin;
1458e6de0934SSatish Balay    hypre_TFree(ptr,HYPRE_MEMORY_HOST);
1459225daaf8SStefano Zampini    PetscFunctionReturn(0);
1460225daaf8SStefano Zampini }
1461