xref: /petsc/src/mat/impls/hypre/mhypre.c (revision c1a070e6be8038c2d0a1a30ecfce97db942ef640)
163c07aadSStefano Zampini 
263c07aadSStefano Zampini /*
363c07aadSStefano Zampini     Creates hypre ijmatrix from PETSc matrix
463c07aadSStefano Zampini */
563c07aadSStefano Zampini #include <petscsys.h>
663c07aadSStefano Zampini #include <petsc/private/matimpl.h>
763c07aadSStefano Zampini #include <../src/mat/impls/hypre/mhypre.h>
863c07aadSStefano Zampini #include <../src/mat/impls/aij/mpi/mpiaij.h>
958968eb6SStefano Zampini #include <../src/vec/vec/impls/hypre/vhyp.h>
1058968eb6SStefano Zampini #include <HYPRE.h>
11*c1a070e6SStefano Zampini #include <HYPRE_utilities.h>
12cd8bc7baSStefano Zampini #include <_hypre_parcsr_ls.h>
1363c07aadSStefano Zampini 
1463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*);
1563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix);
1663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
1763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
1863c07aadSStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool);
1963c07aadSStefano Zampini 
2063c07aadSStefano Zampini #undef __FUNCT__
2163c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate"
2263c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij)
2363c07aadSStefano Zampini {
2463c07aadSStefano Zampini   PetscErrorCode ierr;
2563c07aadSStefano Zampini   PetscInt       i,n_d,n_o;
2663c07aadSStefano Zampini   const PetscInt *ia_d,*ia_o;
2763c07aadSStefano Zampini   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
2863c07aadSStefano Zampini   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
2963c07aadSStefano Zampini 
3063c07aadSStefano Zampini   PetscFunctionBegin;
3163c07aadSStefano Zampini   if (A_d) { /* determine number of nonzero entries in local diagonal part */
3263c07aadSStefano Zampini     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
3363c07aadSStefano Zampini     if (done_d) {
3463c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr);
3563c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
3663c07aadSStefano Zampini         nnz_d[i] = ia_d[i+1] - ia_d[i];
3763c07aadSStefano Zampini       }
3863c07aadSStefano Zampini     }
3963c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
4063c07aadSStefano Zampini   }
4163c07aadSStefano Zampini   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
4263c07aadSStefano Zampini     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
4363c07aadSStefano Zampini     if (done_o) {
4463c07aadSStefano Zampini       ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr);
4563c07aadSStefano Zampini       for (i=0; i<n_o; i++) {
4663c07aadSStefano Zampini         nnz_o[i] = ia_o[i+1] - ia_o[i];
4763c07aadSStefano Zampini       }
4863c07aadSStefano Zampini     }
4963c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
5063c07aadSStefano Zampini   }
5163c07aadSStefano Zampini   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
5263c07aadSStefano Zampini     if (!done_o) { /* only diagonal part */
5363c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
5463c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
5563c07aadSStefano Zampini         nnz_o[i] = 0;
5663c07aadSStefano Zampini       }
5763c07aadSStefano Zampini     }
5863c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o));
5963c07aadSStefano Zampini     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
6063c07aadSStefano Zampini     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
6163c07aadSStefano Zampini   }
6263c07aadSStefano Zampini   PetscFunctionReturn(0);
6363c07aadSStefano Zampini }
6463c07aadSStefano Zampini 
6563c07aadSStefano Zampini #undef __FUNCT__
6663c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_CreateFromMat"
6763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA)
6863c07aadSStefano Zampini {
6963c07aadSStefano Zampini   PetscErrorCode ierr;
7063c07aadSStefano Zampini   PetscInt       rstart,rend,cstart,cend;
7163c07aadSStefano Zampini 
7263c07aadSStefano Zampini   PetscFunctionBegin;
7363c07aadSStefano Zampini   ierr   = MatSetUp(A);CHKERRQ(ierr);
7463c07aadSStefano Zampini   rstart = A->rmap->rstart;
7563c07aadSStefano Zampini   rend   = A->rmap->rend;
7663c07aadSStefano Zampini   cstart = A->cmap->rstart;
7763c07aadSStefano Zampini   cend   = A->cmap->rend;
7863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
7963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
8063c07aadSStefano Zampini   {
8163c07aadSStefano Zampini     PetscBool      same;
8263c07aadSStefano Zampini     Mat            A_d,A_o;
8363c07aadSStefano Zampini     const PetscInt *colmap;
8463c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
8563c07aadSStefano Zampini     if (same) {
8663c07aadSStefano Zampini       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
8763c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
8863c07aadSStefano Zampini       PetscFunctionReturn(0);
8963c07aadSStefano Zampini     }
9063c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
9163c07aadSStefano Zampini     if (same) {
9263c07aadSStefano Zampini       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
9363c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
9463c07aadSStefano Zampini       PetscFunctionReturn(0);
9563c07aadSStefano Zampini     }
9663c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
9763c07aadSStefano Zampini     if (same) {
9863c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
9963c07aadSStefano Zampini       PetscFunctionReturn(0);
10063c07aadSStefano Zampini     }
10163c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
10263c07aadSStefano Zampini     if (same) {
10363c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
10463c07aadSStefano Zampini       PetscFunctionReturn(0);
10563c07aadSStefano Zampini     }
10663c07aadSStefano Zampini   }
10763c07aadSStefano Zampini   PetscFunctionReturn(0);
10863c07aadSStefano Zampini }
10963c07aadSStefano Zampini 
11063c07aadSStefano Zampini #undef __FUNCT__
11163c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixCopy"
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);
13763c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
13863c07aadSStefano Zampini     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
13963c07aadSStefano Zampini   }
14063c07aadSStefano Zampini   PetscFunctionReturn(0);
14163c07aadSStefano Zampini }
14263c07aadSStefano Zampini 
14363c07aadSStefano Zampini #undef __FUNCT__
14463c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ"
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 #undef __FUNCT__
17363c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ"
17463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
17563c07aadSStefano Zampini {
17663c07aadSStefano Zampini   PetscErrorCode        ierr;
17763c07aadSStefano Zampini   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
17863c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag,*poffd;
17963c07aadSStefano Zampini   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
18058968eb6SStefano Zampini   HYPRE_Int             type;
18163c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
18263c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
18363c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag,*hoffd;
18463c07aadSStefano Zampini 
18563c07aadSStefano Zampini   PetscFunctionBegin;
18663c07aadSStefano Zampini   pdiag = (Mat_SeqAIJ*) pA->A->data;
18763c07aadSStefano Zampini   poffd = (Mat_SeqAIJ*) pA->B->data;
18863c07aadSStefano Zampini   /* cstart is only valid for square MPIAIJ layed out in the usual way */
18963c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
19063c07aadSStefano Zampini 
19163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
192ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
193ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
194ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
19563c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
19663c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
19763c07aadSStefano Zampini 
19863c07aadSStefano Zampini   /*
19963c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
20063c07aadSStefano Zampini   */
20163c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
20263c07aadSStefano Zampini   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
20363c07aadSStefano Zampini   jj  = (PetscInt*)hdiag->j;
20463c07aadSStefano Zampini   pjj = (PetscInt*)pdiag->j;
20563c07aadSStefano Zampini   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
20663c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
20763c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
20863c07aadSStefano Zampini   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
20963c07aadSStefano Zampini      If we hacked a hypre a bit more we might be able to avoid this step */
21063c07aadSStefano Zampini   jj  = (PetscInt*) hoffd->j;
21163c07aadSStefano Zampini   pjj = (PetscInt*) poffd->j;
21263c07aadSStefano Zampini   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
21363c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
21463c07aadSStefano Zampini 
215ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
21663c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
21763c07aadSStefano Zampini   PetscFunctionReturn(0);
21863c07aadSStefano Zampini }
21963c07aadSStefano Zampini 
22063c07aadSStefano Zampini #undef __FUNCT__
2212df22349SStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_IS"
2222df22349SStefano Zampini static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B)
2232df22349SStefano Zampini {
2242df22349SStefano Zampini   Mat_HYPRE*             mhA = (Mat_HYPRE*)(A->data);
2252df22349SStefano Zampini   Mat                    lA;
2262df22349SStefano Zampini   ISLocalToGlobalMapping rl2g,cl2g;
2272df22349SStefano Zampini   IS                     is;
2282df22349SStefano Zampini   hypre_ParCSRMatrix     *hA;
2292df22349SStefano Zampini   hypre_CSRMatrix        *hdiag,*hoffd;
2302df22349SStefano Zampini   MPI_Comm               comm;
2312df22349SStefano Zampini   PetscScalar            *hdd,*hod,*aa,*data;
2322df22349SStefano Zampini   PetscInt               *col_map_offd,*hdi,*hdj,*hoi,*hoj;
2332df22349SStefano Zampini   PetscInt               *ii,*jj,*iptr,*jptr;
2342df22349SStefano Zampini   PetscInt               cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N;
23558968eb6SStefano Zampini   HYPRE_Int              type;
2362df22349SStefano Zampini   PetscErrorCode         ierr;
2372df22349SStefano Zampini 
2382df22349SStefano Zampini   PetscFunctionBegin;
239a1787963SStefano Zampini   comm = PetscObjectComm((PetscObject)A);
2402df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type));
2412df22349SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
2422df22349SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA));
2432df22349SStefano Zampini   M     = hypre_ParCSRMatrixGlobalNumRows(hA);
2442df22349SStefano Zampini   N     = hypre_ParCSRMatrixGlobalNumCols(hA);
2452df22349SStefano Zampini   str   = hypre_ParCSRMatrixFirstRowIndex(hA);
2462df22349SStefano Zampini   stc   = hypre_ParCSRMatrixFirstColDiag(hA);
2472df22349SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(hA);
2482df22349SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(hA);
2492df22349SStefano Zampini   dr    = hypre_CSRMatrixNumRows(hdiag);
2502df22349SStefano Zampini   dc    = hypre_CSRMatrixNumCols(hdiag);
2512df22349SStefano Zampini   nnz   = hypre_CSRMatrixNumNonzeros(hdiag);
2522df22349SStefano Zampini   hdi   = hypre_CSRMatrixI(hdiag);
2532df22349SStefano Zampini   hdj   = hypre_CSRMatrixJ(hdiag);
2542df22349SStefano Zampini   hdd   = hypre_CSRMatrixData(hdiag);
2552df22349SStefano Zampini   oc    = hypre_CSRMatrixNumCols(hoffd);
2562df22349SStefano Zampini   nnz  += hypre_CSRMatrixNumNonzeros(hoffd);
2572df22349SStefano Zampini   hoi   = hypre_CSRMatrixI(hoffd);
2582df22349SStefano Zampini   hoj   = hypre_CSRMatrixJ(hoffd);
2592df22349SStefano Zampini   hod   = hypre_CSRMatrixData(hoffd);
2602df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
2612df22349SStefano Zampini     PetscInt *aux;
2622df22349SStefano Zampini 
2632df22349SStefano Zampini     /* generate l2g maps for rows and cols */
2642df22349SStefano Zampini     ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr);
2652df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
2662df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2672df22349SStefano Zampini     col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
2682df22349SStefano Zampini     ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr);
2692df22349SStefano Zampini     for (i=0; i<dc; i++) aux[i] = i+stc;
2702df22349SStefano Zampini     for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i];
2712df22349SStefano Zampini     ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
2722df22349SStefano Zampini     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
2732df22349SStefano Zampini     ierr = ISDestroy(&is);CHKERRQ(ierr);
2742df22349SStefano Zampini     /* create MATIS object */
2752df22349SStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
2762df22349SStefano Zampini     ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr);
2772df22349SStefano Zampini     ierr = MatSetType(*B,MATIS);CHKERRQ(ierr);
2782df22349SStefano Zampini     ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr);
2792df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
2802df22349SStefano Zampini     ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
2812df22349SStefano Zampini 
2822df22349SStefano Zampini     /* allocate CSR for local matrix */
2832df22349SStefano Zampini     ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr);
2842df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr);
2852df22349SStefano Zampini     ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr);
2862df22349SStefano Zampini   } else {
2872df22349SStefano Zampini     PetscInt  nr;
2882df22349SStefano Zampini     PetscBool done;
2892df22349SStefano Zampini     ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr);
2902df22349SStefano Zampini     ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr);
2912df22349SStefano 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);
2922df22349SStefano 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);
2932df22349SStefano Zampini     ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr);
2942df22349SStefano Zampini   }
2952df22349SStefano Zampini   /* merge local matrices */
2962df22349SStefano Zampini   ii   = iptr;
2972df22349SStefano Zampini   jj   = jptr;
2982df22349SStefano Zampini   aa   = data;
2992df22349SStefano Zampini   *ii  = *(hdi++) + *(hoi++);
3002df22349SStefano Zampini   for (jd=0,jo=0,cum=0; *ii<nnz; cum++) {
3012df22349SStefano Zampini     PetscScalar *aold = aa;
3022df22349SStefano Zampini     PetscInt    *jold = jj,nc = jd+jo;
3032df22349SStefano Zampini     for (; jd<*hdi; jd++) { *jj++ = *hdj++;      *aa++ = *hdd++; }
3042df22349SStefano Zampini     for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
3052df22349SStefano Zampini     *(++ii) = *(hdi++) + *(hoi++);
3062df22349SStefano Zampini     ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr);
3072df22349SStefano Zampini   }
3082df22349SStefano Zampini   for (; cum<dr; cum++) *(++ii) = nnz;
3092df22349SStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
310a033916dSStefano Zampini     Mat_SeqAIJ* a;
311a033916dSStefano Zampini 
3122df22349SStefano Zampini     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr);
3132df22349SStefano Zampini     ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr);
314a033916dSStefano Zampini     /* hack SeqAIJ */
315a033916dSStefano Zampini     a          = (Mat_SeqAIJ*)(lA->data);
316a033916dSStefano Zampini     a->free_a  = PETSC_TRUE;
317a033916dSStefano Zampini     a->free_ij = PETSC_TRUE;
3182df22349SStefano Zampini     ierr = MatDestroy(&lA);CHKERRQ(ierr);
3192df22349SStefano Zampini   }
3202df22349SStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3212df22349SStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
3222df22349SStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
3232df22349SStefano Zampini     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
3242df22349SStefano Zampini   }
3252df22349SStefano Zampini   PetscFunctionReturn(0);
3262df22349SStefano Zampini }
3272df22349SStefano Zampini 
3282df22349SStefano Zampini #undef __FUNCT__
32963c07aadSStefano Zampini #define __FUNCT__ "MatConvert_AIJ_HYPRE"
33063c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
33163c07aadSStefano Zampini {
33263c07aadSStefano Zampini   Mat_HYPRE      *hB;
33363c07aadSStefano Zampini   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
33463c07aadSStefano Zampini   PetscErrorCode ierr;
33563c07aadSStefano Zampini 
33663c07aadSStefano Zampini   PetscFunctionBegin;
33763c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
33863c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
33963c07aadSStefano Zampini     /* always destroy the old matrix and create a new memory;
34063c07aadSStefano Zampini        hope this does not churn the memory too much. The problem
34163c07aadSStefano Zampini        is I do not know if it is possible to put the matrix back to
34263c07aadSStefano Zampini        its initial state so that we can directly copy the values
34363c07aadSStefano Zampini        the second time through. */
34463c07aadSStefano Zampini     hB = (Mat_HYPRE*)((*B)->data);
34563c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
34663c07aadSStefano Zampini   } else {
34763c07aadSStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
34863c07aadSStefano Zampini     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
34963c07aadSStefano Zampini     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
35063c07aadSStefano Zampini     ierr = MatSetUp(*B);CHKERRQ(ierr);
35163c07aadSStefano Zampini     hB   = (Mat_HYPRE*)((*B)->data);
35263c07aadSStefano Zampini   }
35363c07aadSStefano Zampini   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
35463c07aadSStefano Zampini   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
35563c07aadSStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
35663c07aadSStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
35763c07aadSStefano Zampini   PetscFunctionReturn(0);
35863c07aadSStefano Zampini }
35963c07aadSStefano Zampini 
36063c07aadSStefano Zampini #undef __FUNCT__
36163c07aadSStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_AIJ"
362ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
36363c07aadSStefano Zampini {
36463c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
36563c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
36663c07aadSStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
36763c07aadSStefano Zampini   MPI_Comm           comm;
36863c07aadSStefano Zampini   PetscScalar        *da,*oa,*aptr;
36963c07aadSStefano Zampini   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
37063c07aadSStefano Zampini   PetscInt           i,dnnz,onnz,m,n;
37158968eb6SStefano Zampini   HYPRE_Int          type;
37263c07aadSStefano Zampini   PetscMPIInt        size;
37363c07aadSStefano Zampini   PetscErrorCode     ierr;
37463c07aadSStefano Zampini 
37563c07aadSStefano Zampini   PetscFunctionBegin;
37663c07aadSStefano Zampini   comm = PetscObjectComm((PetscObject)A);
37763c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
37863c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
37963c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
38063c07aadSStefano Zampini     PetscBool ismpiaij,isseqaij;
38163c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
38263c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
38363c07aadSStefano Zampini     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
38463c07aadSStefano Zampini   }
38563c07aadSStefano Zampini   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
38663c07aadSStefano Zampini 
38763c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
38863c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(parcsr);
38963c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(parcsr);
39063c07aadSStefano Zampini   m     = hypre_CSRMatrixNumRows(hdiag);
39163c07aadSStefano Zampini   n     = hypre_CSRMatrixNumCols(hdiag);
39263c07aadSStefano Zampini   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
39363c07aadSStefano Zampini   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
39463c07aadSStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
39563c07aadSStefano Zampini     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
39663c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
39763c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
39863c07aadSStefano Zampini   } else {
39963c07aadSStefano Zampini     PetscInt  nr;
40063c07aadSStefano Zampini     PetscBool done;
40163c07aadSStefano Zampini     if (size > 1) {
40263c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
40363c07aadSStefano Zampini 
40463c07aadSStefano Zampini       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
40563c07aadSStefano 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);
40663c07aadSStefano 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);
40763c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
40863c07aadSStefano Zampini     } else {
40963c07aadSStefano Zampini       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
41063c07aadSStefano Zampini       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
41163c07aadSStefano 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);
41263c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
41363c07aadSStefano Zampini     }
41463c07aadSStefano Zampini   }
41563c07aadSStefano Zampini   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
41663c07aadSStefano Zampini   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
41763c07aadSStefano Zampini   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
41863c07aadSStefano Zampini   iptr = djj;
41963c07aadSStefano Zampini   aptr = da;
42063c07aadSStefano Zampini   for (i=0; i<m; i++) {
42163c07aadSStefano Zampini     PetscInt nc = dii[i+1]-dii[i];
42263c07aadSStefano Zampini     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
42363c07aadSStefano Zampini     iptr += nc;
42463c07aadSStefano Zampini     aptr += nc;
42563c07aadSStefano Zampini   }
42663c07aadSStefano Zampini   if (size > 1) {
42763c07aadSStefano Zampini     PetscInt *offdj,*coffd;
42863c07aadSStefano Zampini 
42963c07aadSStefano Zampini     if (reuse != MAT_REUSE_MATRIX) {
43063c07aadSStefano Zampini       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
43163c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
43263c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
43363c07aadSStefano Zampini     } else {
43463c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
43563c07aadSStefano Zampini       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
43663c07aadSStefano Zampini       PetscBool  done;
43763c07aadSStefano Zampini 
43863c07aadSStefano Zampini       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
43963c07aadSStefano 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);
44063c07aadSStefano 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);
44163c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
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     }
45663c07aadSStefano Zampini     if (reuse != MAT_REUSE_MATRIX) {
45763c07aadSStefano Zampini       Mat_MPIAIJ *b;
45863c07aadSStefano Zampini       Mat_SeqAIJ *d,*o;
45963c07aadSStefano Zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,
46063c07aadSStefano Zampini                                             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;
46963c07aadSStefano Zampini     }
47063c07aadSStefano Zampini   } else if (reuse != MAT_REUSE_MATRIX) {
47163c07aadSStefano Zampini     Mat_SeqAIJ* b;
47263c07aadSStefano Zampini     ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
47363c07aadSStefano Zampini     /* hack SeqAIJ */
47463c07aadSStefano Zampini     b          = (Mat_SeqAIJ*)((*B)->data);
47563c07aadSStefano Zampini     b->free_a  = PETSC_TRUE;
47663c07aadSStefano Zampini     b->free_ij = PETSC_TRUE;
47763c07aadSStefano Zampini   }
47863c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
47963c07aadSStefano Zampini     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
48063c07aadSStefano Zampini   }
48163c07aadSStefano Zampini   PetscFunctionReturn(0);
48263c07aadSStefano Zampini }
48363c07aadSStefano Zampini 
48463c07aadSStefano Zampini #undef __FUNCT__
485*c1a070e6SStefano Zampini #define __FUNCT__ "MatPtAP_AIJ_HYPRE"
486*c1a070e6SStefano Zampini static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
487*c1a070e6SStefano Zampini {
488*c1a070e6SStefano Zampini   Mat_HYPRE          *hP = (Mat_HYPRE*)P->data;
489*c1a070e6SStefano Zampini   hypre_ParCSRMatrix *tA,*Pparcsr,*ptapparcsr;
490*c1a070e6SStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
491*c1a070e6SStefano Zampini   Mat_SeqAIJ         *diag,*offd;
492*c1a070e6SStefano Zampini   PetscInt           *garray,noffd,dnnz,onnz,*row_starts,*col_starts;
493*c1a070e6SStefano Zampini   HYPRE_Int          type,P_owns_col_starts;
494*c1a070e6SStefano Zampini   PetscBool          ismpiaij,isseqaij;
495*c1a070e6SStefano Zampini   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
496*c1a070e6SStefano Zampini   PetscErrorCode     ierr;
497*c1a070e6SStefano Zampini 
498*c1a070e6SStefano Zampini   PetscFunctionBegin;
499*c1a070e6SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
500*c1a070e6SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
501*c1a070e6SStefano Zampini   if (scall == MAT_REUSE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_REUSE_MATRIX");
502*c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
503*c1a070e6SStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
504*c1a070e6SStefano Zampini   if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type);
505*c1a070e6SStefano Zampini   if (ismpiaij) {
506*c1a070e6SStefano Zampini     Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data);
507*c1a070e6SStefano Zampini 
508*c1a070e6SStefano Zampini     diag   = (Mat_SeqAIJ*)a->A->data;
509*c1a070e6SStefano Zampini     offd   = (Mat_SeqAIJ*)a->B->data;
510*c1a070e6SStefano Zampini     garray = a->garray;
511*c1a070e6SStefano Zampini     noffd  = a->B->cmap->N;
512*c1a070e6SStefano Zampini     dnnz   = diag->nz;
513*c1a070e6SStefano Zampini     onnz   = offd->nz;
514*c1a070e6SStefano Zampini   } else {
515*c1a070e6SStefano Zampini     diag    = (Mat_SeqAIJ*)A->data;
516*c1a070e6SStefano Zampini     offd    = NULL;
517*c1a070e6SStefano Zampini     garray  = NULL;
518*c1a070e6SStefano Zampini     noffd   = 0;
519*c1a070e6SStefano Zampini     dnnz    = diag->nz;
520*c1a070e6SStefano Zampini     onnz    = 0;
521*c1a070e6SStefano Zampini   }
522*c1a070e6SStefano Zampini   /* create a temporary ParCSR */
523*c1a070e6SStefano Zampini   if (HYPRE_AssumedPartitionCheck()) {
524*c1a070e6SStefano Zampini    PetscMPIInt myid;
525*c1a070e6SStefano Zampini 
526*c1a070e6SStefano Zampini    ierr       = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
527*c1a070e6SStefano Zampini    row_starts = A->rmap->range + myid;
528*c1a070e6SStefano Zampini    col_starts = A->cmap->range + myid;
529*c1a070e6SStefano Zampini   } else {
530*c1a070e6SStefano Zampini    row_starts = A->rmap->range;
531*c1a070e6SStefano Zampini    col_starts = A->cmap->range;
532*c1a070e6SStefano Zampini   }
533*c1a070e6SStefano Zampini   tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,row_starts,col_starts,noffd,dnnz,onnz);
534*c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(tA,0);
535*c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(tA,0);
536*c1a070e6SStefano Zampini 
537*c1a070e6SStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(tA);
538*c1a070e6SStefano Zampini   hypre_CSRMatrixI(hdiag)           = diag->i;
539*c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hdiag)           = diag->j;
540*c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)        = diag->a;
541*c1a070e6SStefano Zampini   hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz;
542*c1a070e6SStefano Zampini   hypre_CSRMatrixSetRownnz(hdiag);
543*c1a070e6SStefano Zampini   hypre_CSRMatrixSetDataOwner(hdiag,0);
544*c1a070e6SStefano Zampini 
545*c1a070e6SStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(tA);
546*c1a070e6SStefano Zampini   if (offd) {
547*c1a070e6SStefano Zampini     hypre_CSRMatrixI(hoffd)           = offd->i;
548*c1a070e6SStefano Zampini     hypre_CSRMatrixJ(hoffd)           = offd->j;
549*c1a070e6SStefano Zampini     hypre_CSRMatrixData(hoffd)        = offd->a;
550*c1a070e6SStefano Zampini     hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz;
551*c1a070e6SStefano Zampini     hypre_CSRMatrixSetRownnz(hoffd);
552*c1a070e6SStefano Zampini     hypre_CSRMatrixSetDataOwner(hoffd,0);
553*c1a070e6SStefano Zampini     hypre_ParCSRMatrixSetNumNonzeros(tA);
554*c1a070e6SStefano Zampini     hypre_ParCSRMatrixColMapOffd(tA) = garray;
555*c1a070e6SStefano Zampini   }
556*c1a070e6SStefano Zampini 
557*c1a070e6SStefano Zampini   /* call RAP from BoomerAMG */
558*c1a070e6SStefano Zampini   /* Warning: hypre_BoomerAMGBuildCoarseOperator steals the col_starts
559*c1a070e6SStefano Zampini      from Pparcsr (even if it does not own them)! */
560*c1a070e6SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
561*c1a070e6SStefano Zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(Pparcsr);
562*c1a070e6SStefano Zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
563*c1a070e6SStefano Zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(Pparcsr,tA,Pparcsr,&ptapparcsr));
564*c1a070e6SStefano Zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
565*c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(ptapparcsr,0);
566*c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(ptapparcsr,0);
567*c1a070e6SStefano Zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(Pparcsr, 1);
568*c1a070e6SStefano Zampini 
569*c1a070e6SStefano Zampini   /* create MatHYPRE */
570*c1a070e6SStefano Zampini   ierr = MatHYPRECreateFromParCSR(ptapparcsr,PETSC_OWN_POINTER,C);CHKERRQ(ierr);
571*c1a070e6SStefano Zampini 
572*c1a070e6SStefano Zampini   /* set pointers to NULL before destroying tA */
573*c1a070e6SStefano Zampini   hypre_CSRMatrixI(hdiag)          = NULL;
574*c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hdiag)          = NULL;
575*c1a070e6SStefano Zampini   hypre_CSRMatrixData(hdiag)       = NULL;
576*c1a070e6SStefano Zampini   hypre_CSRMatrixI(hoffd)          = NULL;
577*c1a070e6SStefano Zampini   hypre_CSRMatrixJ(hoffd)          = NULL;
578*c1a070e6SStefano Zampini   hypre_CSRMatrixData(hoffd)       = NULL;
579*c1a070e6SStefano Zampini   hypre_ParCSRMatrixColMapOffd(tA) = NULL;
580*c1a070e6SStefano Zampini   hypre_ParCSRMatrixDestroy(tA);
581*c1a070e6SStefano Zampini   PetscFunctionReturn(0);
582*c1a070e6SStefano Zampini }
583*c1a070e6SStefano Zampini 
584*c1a070e6SStefano Zampini #undef __FUNCT__
585cd8bc7baSStefano Zampini #define __FUNCT__ "MatPtAP_HYPRE_HYPRE"
586cd8bc7baSStefano Zampini static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
587cd8bc7baSStefano Zampini {
588cd8bc7baSStefano Zampini   hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr;
589cd8bc7baSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data, *hP = (Mat_HYPRE*)P->data;
590*c1a070e6SStefano Zampini   HYPRE_Int          type,P_owns_col_starts;
591cd8bc7baSStefano Zampini   PetscErrorCode     ierr;
592cd8bc7baSStefano Zampini 
593cd8bc7baSStefano Zampini   PetscFunctionBegin;
594cd8bc7baSStefano Zampini   if (scall == MAT_REUSE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Unsupported MAT_REUSE_MATRIX");
595cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
596cd8bc7baSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
597cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
598cd8bc7baSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
599cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
600cd8bc7baSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
601*c1a070e6SStefano Zampini 
602*c1a070e6SStefano Zampini   /* call RAP from BoomerAMG */
603*c1a070e6SStefano Zampini   /* Warning: hypre_BoomerAMGBuildCoarseOperator steals the col_starts
604*c1a070e6SStefano Zampini      from Pparcsr (even if it does not own them)! */
605*c1a070e6SStefano Zampini   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(Pparcsr);
606*c1a070e6SStefano Zampini   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
607cd8bc7baSStefano Zampini   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr));
608cd8bc7baSStefano Zampini   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
609*c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetRowStartsOwner(ptapparcsr,0);
610*c1a070e6SStefano Zampini   hypre_ParCSRMatrixSetColStartsOwner(ptapparcsr,0);
611*c1a070e6SStefano Zampini   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(Pparcsr, 1);
612*c1a070e6SStefano Zampini 
613*c1a070e6SStefano Zampini   /* create MatHYPRE */
614*c1a070e6SStefano Zampini   ierr = MatHYPRECreateFromParCSR(ptapparcsr,PETSC_OWN_POINTER,C);CHKERRQ(ierr);
615cd8bc7baSStefano Zampini   PetscFunctionReturn(0);
616cd8bc7baSStefano Zampini }
617cd8bc7baSStefano Zampini 
618cd8bc7baSStefano Zampini #undef __FUNCT__
61963c07aadSStefano Zampini #define __FUNCT__ "MatMultTranspose_HYPRE"
620ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
62163c07aadSStefano Zampini {
62263c07aadSStefano Zampini   PetscErrorCode ierr;
62363c07aadSStefano Zampini 
62463c07aadSStefano Zampini   PetscFunctionBegin;
62563c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
62663c07aadSStefano Zampini   PetscFunctionReturn(0);
62763c07aadSStefano Zampini }
62863c07aadSStefano Zampini 
62963c07aadSStefano Zampini #undef __FUNCT__
63063c07aadSStefano Zampini #define __FUNCT__ "MatMult_HYPRE"
631ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
63263c07aadSStefano Zampini {
63363c07aadSStefano Zampini   PetscErrorCode ierr;
63463c07aadSStefano Zampini 
63563c07aadSStefano Zampini   PetscFunctionBegin;
63663c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
63763c07aadSStefano Zampini   PetscFunctionReturn(0);
63863c07aadSStefano Zampini }
63963c07aadSStefano Zampini 
64063c07aadSStefano Zampini #undef __FUNCT__
64163c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_MultKernel_Private"
642ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
64363c07aadSStefano Zampini {
64463c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
64563c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
64663c07aadSStefano Zampini   hypre_ParVector    *hx,*hy;
64763c07aadSStefano Zampini   PetscScalar        *ax,*ay,*sax,*say;
64863c07aadSStefano Zampini   PetscErrorCode     ierr;
64963c07aadSStefano Zampini 
65063c07aadSStefano Zampini   PetscFunctionBegin;
65163c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
65263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
65363c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
65463c07aadSStefano Zampini   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
65563c07aadSStefano Zampini   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
65663c07aadSStefano Zampini   if (trans) {
65758968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ay,say);
65858968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax);
65963c07aadSStefano Zampini     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
66058968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,say,ay);
66158968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax);
66263c07aadSStefano Zampini   } else {
66358968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax);
66458968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,ay,say);
66563c07aadSStefano Zampini     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
66658968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax);
66758968eb6SStefano Zampini     VecHYPRE_ParVectorReplacePointer(hA->b,say,ay);
66863c07aadSStefano Zampini   }
66963c07aadSStefano Zampini   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
67063c07aadSStefano Zampini   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
67163c07aadSStefano Zampini   PetscFunctionReturn(0);
67263c07aadSStefano Zampini }
67363c07aadSStefano Zampini 
67463c07aadSStefano Zampini #undef __FUNCT__
67563c07aadSStefano Zampini #define __FUNCT__ "MatDestroy_HYPRE"
676ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A)
67763c07aadSStefano Zampini {
67863c07aadSStefano Zampini   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
67963c07aadSStefano Zampini   PetscErrorCode ierr;
68063c07aadSStefano Zampini 
68163c07aadSStefano Zampini   PetscFunctionBegin;
68263c07aadSStefano Zampini   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
68363c07aadSStefano Zampini   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
684978814f1SStefano Zampini   if (hA->ij) {
685978814f1SStefano Zampini     if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL;
686978814f1SStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
687978814f1SStefano Zampini   }
68863c07aadSStefano Zampini   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); }
68963c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
6902df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr);
691*c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr);
692*c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr);
69363c07aadSStefano Zampini   ierr = PetscFree(A->data);CHKERRQ(ierr);
69463c07aadSStefano Zampini   PetscFunctionReturn(0);
69563c07aadSStefano Zampini }
69663c07aadSStefano Zampini 
69763c07aadSStefano Zampini #undef __FUNCT__
69863c07aadSStefano Zampini #define __FUNCT__ "MatSetUp_HYPRE"
699ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A)
70063c07aadSStefano Zampini {
70163c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
70263c07aadSStefano Zampini   Vec                x,b;
70363c07aadSStefano Zampini   PetscErrorCode     ierr;
70463c07aadSStefano Zampini 
70563c07aadSStefano Zampini   PetscFunctionBegin;
70663c07aadSStefano Zampini   if (hA->x) PetscFunctionReturn(0);
70763c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
70863c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
70963c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
71063c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
71163c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
71263c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
71363c07aadSStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
71463c07aadSStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
71563c07aadSStefano Zampini   PetscFunctionReturn(0);
71663c07aadSStefano Zampini }
71763c07aadSStefano Zampini 
71863c07aadSStefano Zampini #undef __FUNCT__
71963c07aadSStefano Zampini #define __FUNCT__ "MatAssemblyEnd_HYPRE"
720ea9daf28SStefano Zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
72163c07aadSStefano Zampini {
72263c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
72363c07aadSStefano Zampini   PetscErrorCode     ierr;
72463c07aadSStefano Zampini 
72563c07aadSStefano Zampini   PetscFunctionBegin;
72663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
72763c07aadSStefano Zampini   PetscFunctionReturn(0);
72863c07aadSStefano Zampini }
72963c07aadSStefano Zampini 
730978814f1SStefano Zampini #undef __FUNCT__
731978814f1SStefano Zampini #define __FUNCT__ "MatHYPRECreateFromParCSR"
732978814f1SStefano Zampini PETSC_EXTERN PetscErrorCode MatHYPRECreateFromParCSR(void *vparcsr, PetscCopyMode copymode, Mat* A)
733978814f1SStefano Zampini {
734978814f1SStefano Zampini   Mat_HYPRE             *hA;
73563c07aadSStefano Zampini   hypre_ParCSRMatrix    *parcsr;
736978814f1SStefano Zampini   MPI_Comm              comm;
737978814f1SStefano Zampini   PetscInt              rstart,rend,cstart,cend,M,N;
738978814f1SStefano Zampini   PetscErrorCode        ierr;
739978814f1SStefano Zampini 
740978814f1SStefano Zampini   PetscFunctionBegin;
741978814f1SStefano Zampini   parcsr = (hypre_ParCSRMatrix *)vparcsr;
742978814f1SStefano Zampini   comm   = hypre_ParCSRMatrixComm(parcsr);
743978814f1SStefano Zampini   if (copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES");
744978814f1SStefano Zampini 
745978814f1SStefano Zampini   /* access ParCSRMatrix */
746978814f1SStefano Zampini   rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr);
747978814f1SStefano Zampini   rend   = hypre_ParCSRMatrixLastRowIndex(parcsr);
748978814f1SStefano Zampini   cstart = hypre_ParCSRMatrixFirstColDiag(parcsr);
749978814f1SStefano Zampini   cend   = hypre_ParCSRMatrixLastColDiag(parcsr);
750978814f1SStefano Zampini   M      = hypre_ParCSRMatrixGlobalNumRows(parcsr);
751978814f1SStefano Zampini   N      = hypre_ParCSRMatrixGlobalNumCols(parcsr);
752978814f1SStefano Zampini 
753978814f1SStefano Zampini   /* create PETSc matrix with MatHYPRE */
754978814f1SStefano Zampini   ierr = MatCreate(comm,A);CHKERRQ(ierr);
755978814f1SStefano Zampini   ierr = MatSetSizes(*A,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr);
756978814f1SStefano Zampini   ierr = MatSetType(*A,MATHYPRE);CHKERRQ(ierr);
757978814f1SStefano Zampini   ierr = MatSetUp(*A);CHKERRQ(ierr);
758978814f1SStefano Zampini   hA   = (Mat_HYPRE*)((*A)->data);
759978814f1SStefano Zampini 
760978814f1SStefano Zampini   /* create HYPRE_IJMatrix */
761978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
762978814f1SStefano Zampini 
763978814f1SStefano Zampini   /* set ParCSR object */
764978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
765978814f1SStefano Zampini   hypre_IJMatrixObject(hA->ij) = parcsr;
766978814f1SStefano Zampini 
767978814f1SStefano Zampini   /* set assembled flag */
768978814f1SStefano Zampini   hypre_IJMatrixAssembleFlag(hA->ij) = 1;
769978814f1SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
770978814f1SStefano Zampini 
771978814f1SStefano Zampini   /* prevent from freeing the pointer */
772978814f1SStefano Zampini   if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE;
773978814f1SStefano Zampini 
774978814f1SStefano Zampini   ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
775978814f1SStefano Zampini   ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
776978814f1SStefano Zampini   PetscFunctionReturn(0);
777978814f1SStefano Zampini }
778978814f1SStefano Zampini 
77963c07aadSStefano Zampini #undef __FUNCT__
78063c07aadSStefano Zampini #define __FUNCT__ "MatCreate_HYPRE"
78163c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
78263c07aadSStefano Zampini {
78363c07aadSStefano Zampini   Mat_HYPRE      *hB;
78463c07aadSStefano Zampini   PetscErrorCode ierr;
78563c07aadSStefano Zampini 
78663c07aadSStefano Zampini   PetscFunctionBegin;
78763c07aadSStefano Zampini   ierr           = PetscNewLog(B,&hB);CHKERRQ(ierr);
788978814f1SStefano Zampini   hB->inner_free = PETSC_TRUE;
789978814f1SStefano Zampini 
79063c07aadSStefano Zampini   B->data       = (void*)hB;
79163c07aadSStefano Zampini   B->rmap->bs   = 1;
79263c07aadSStefano Zampini   B->assembled  = PETSC_FALSE;
79363c07aadSStefano Zampini 
79463c07aadSStefano Zampini   B->ops->mult          = MatMult_HYPRE;
79563c07aadSStefano Zampini   B->ops->multtranspose = MatMultTranspose_HYPRE;
79663c07aadSStefano Zampini   B->ops->setup         = MatSetUp_HYPRE;
79763c07aadSStefano Zampini   B->ops->destroy       = MatDestroy_HYPRE;
79863c07aadSStefano Zampini   B->ops->assemblyend   = MatAssemblyEnd_HYPRE;
799cd8bc7baSStefano Zampini   B->ops->ptap          = MatPtAP_HYPRE_HYPRE;
80063c07aadSStefano Zampini 
80163c07aadSStefano Zampini   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
80263c07aadSStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
80363c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
8042df22349SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr);
805*c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
806*c1a070e6SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
80763c07aadSStefano Zampini   PetscFunctionReturn(0);
80863c07aadSStefano Zampini }
80963c07aadSStefano Zampini 
81063c07aadSStefano Zampini #if 0
81163c07aadSStefano Zampini /*
81263c07aadSStefano Zampini     Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format
81363c07aadSStefano Zampini 
81463c07aadSStefano Zampini     This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first
81563c07aadSStefano Zampini     which will corrupt the PETSc data structure if we did this. Need a work around to this problem.
81663c07aadSStefano Zampini */
81763c07aadSStefano Zampini #include <_hypre_IJ_mv.h>
81863c07aadSStefano Zampini #include <HYPRE_IJ_mv.h>
81963c07aadSStefano Zampini 
82063c07aadSStefano Zampini #undef __FUNCT__
82163c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixLink"
82263c07aadSStefano Zampini PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij)
82363c07aadSStefano Zampini {
82463c07aadSStefano Zampini   PetscErrorCode        ierr;
82563c07aadSStefano Zampini   PetscInt              rstart,rend,cstart,cend;
82663c07aadSStefano Zampini   PetscBool             flg;
82763c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
82863c07aadSStefano Zampini 
82963c07aadSStefano Zampini   PetscFunctionBegin;
83063c07aadSStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
83163c07aadSStefano Zampini   PetscValidType(A,1);
83263c07aadSStefano Zampini   PetscValidPointer(ij,2);
83363c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
83463c07aadSStefano Zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices");
83563c07aadSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
83663c07aadSStefano Zampini 
83763c07aadSStefano Zampini   ierr   = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
83863c07aadSStefano Zampini   rstart = A->rmap->rstart;
83963c07aadSStefano Zampini   rend   = A->rmap->rend;
84063c07aadSStefano Zampini   cstart = A->cmap->rstart;
84163c07aadSStefano Zampini   cend   = A->cmap->rend;
84263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij));
84363c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
84463c07aadSStefano Zampini 
84563c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij));
84663c07aadSStefano Zampini   PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij));
84763c07aadSStefano Zampini 
84863c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
84963c07aadSStefano Zampini 
85063c07aadSStefano Zampini   /* this is the Hack part where we monkey directly with the hypre datastructures */
85163c07aadSStefano Zampini 
85263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij));
85363c07aadSStefano Zampini   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
85463c07aadSStefano Zampini   PetscFunctionReturn(0);
85563c07aadSStefano Zampini }
85663c07aadSStefano Zampini #endif
857