xref: /petsc/src/mat/impls/hypre/mhypre.c (revision ea9daf28f8d450bf21161ac4a6ee3117b9487a48)
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>
963c07aadSStefano Zampini #include <HYPRE_struct_mv.h>
1063c07aadSStefano Zampini #include <HYPRE_struct_ls.h>
1163c07aadSStefano Zampini #include <_hypre_struct_mv.h>
1263c07aadSStefano Zampini 
1363c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*);
1463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix);
1563c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
1663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
1763c07aadSStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool);
1863c07aadSStefano Zampini 
1963c07aadSStefano Zampini #undef __FUNCT__
2063c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate"
2163c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij)
2263c07aadSStefano Zampini {
2363c07aadSStefano Zampini   PetscErrorCode ierr;
2463c07aadSStefano Zampini   PetscInt       i,n_d,n_o;
2563c07aadSStefano Zampini   const PetscInt *ia_d,*ia_o;
2663c07aadSStefano Zampini   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
2763c07aadSStefano Zampini   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
2863c07aadSStefano Zampini 
2963c07aadSStefano Zampini   PetscFunctionBegin;
3063c07aadSStefano Zampini   if (A_d) { /* determine number of nonzero entries in local diagonal part */
3163c07aadSStefano Zampini     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
3263c07aadSStefano Zampini     if (done_d) {
3363c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr);
3463c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
3563c07aadSStefano Zampini         nnz_d[i] = ia_d[i+1] - ia_d[i];
3663c07aadSStefano Zampini       }
3763c07aadSStefano Zampini     }
3863c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
3963c07aadSStefano Zampini   }
4063c07aadSStefano Zampini   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
4163c07aadSStefano Zampini     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
4263c07aadSStefano Zampini     if (done_o) {
4363c07aadSStefano Zampini       ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr);
4463c07aadSStefano Zampini       for (i=0; i<n_o; i++) {
4563c07aadSStefano Zampini         nnz_o[i] = ia_o[i+1] - ia_o[i];
4663c07aadSStefano Zampini       }
4763c07aadSStefano Zampini     }
4863c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
4963c07aadSStefano Zampini   }
5063c07aadSStefano Zampini   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
5163c07aadSStefano Zampini     if (!done_o) { /* only diagonal part */
5263c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
5363c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
5463c07aadSStefano Zampini         nnz_o[i] = 0;
5563c07aadSStefano Zampini       }
5663c07aadSStefano Zampini     }
5763c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o));
5863c07aadSStefano Zampini     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
5963c07aadSStefano Zampini     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
6063c07aadSStefano Zampini   }
6163c07aadSStefano Zampini   PetscFunctionReturn(0);
6263c07aadSStefano Zampini }
6363c07aadSStefano Zampini 
6463c07aadSStefano Zampini #undef __FUNCT__
6563c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_CreateFromMat"
6663c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA)
6763c07aadSStefano Zampini {
6863c07aadSStefano Zampini   PetscErrorCode ierr;
6963c07aadSStefano Zampini   PetscInt       rstart,rend,cstart,cend;
7063c07aadSStefano Zampini 
7163c07aadSStefano Zampini   PetscFunctionBegin;
7263c07aadSStefano Zampini   ierr   = MatSetUp(A);CHKERRQ(ierr);
7363c07aadSStefano Zampini   rstart = A->rmap->rstart;
7463c07aadSStefano Zampini   rend   = A->rmap->rend;
7563c07aadSStefano Zampini   cstart = A->cmap->rstart;
7663c07aadSStefano Zampini   cend   = A->cmap->rend;
7763c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
7863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
7963c07aadSStefano Zampini   {
8063c07aadSStefano Zampini     PetscBool      same;
8163c07aadSStefano Zampini     Mat            A_d,A_o;
8263c07aadSStefano Zampini     const PetscInt *colmap;
8363c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
8463c07aadSStefano Zampini     if (same) {
8563c07aadSStefano Zampini       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
8663c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
8763c07aadSStefano Zampini       PetscFunctionReturn(0);
8863c07aadSStefano Zampini     }
8963c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
9063c07aadSStefano Zampini     if (same) {
9163c07aadSStefano Zampini       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
9263c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
9363c07aadSStefano Zampini       PetscFunctionReturn(0);
9463c07aadSStefano Zampini     }
9563c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
9663c07aadSStefano Zampini     if (same) {
9763c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
9863c07aadSStefano Zampini       PetscFunctionReturn(0);
9963c07aadSStefano Zampini     }
10063c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
10163c07aadSStefano Zampini     if (same) {
10263c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
10363c07aadSStefano Zampini       PetscFunctionReturn(0);
10463c07aadSStefano Zampini     }
10563c07aadSStefano Zampini   }
10663c07aadSStefano Zampini   PetscFunctionReturn(0);
10763c07aadSStefano Zampini }
10863c07aadSStefano Zampini 
10963c07aadSStefano Zampini #undef __FUNCT__
11063c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixCopy"
11163c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij)
11263c07aadSStefano Zampini {
11363c07aadSStefano Zampini   PetscErrorCode    ierr;
11463c07aadSStefano Zampini   PetscInt          i,rstart,rend,ncols,nr,nc;
11563c07aadSStefano Zampini   const PetscScalar *values;
11663c07aadSStefano Zampini   const PetscInt    *cols;
11763c07aadSStefano Zampini   PetscBool         flg;
11863c07aadSStefano Zampini 
11963c07aadSStefano Zampini   PetscFunctionBegin;
12063c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
12163c07aadSStefano Zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
12263c07aadSStefano Zampini   if (flg && nr == nc) {
12363c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
12463c07aadSStefano Zampini     PetscFunctionReturn(0);
12563c07aadSStefano Zampini   }
12663c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
12763c07aadSStefano Zampini   if (flg) {
12863c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
12963c07aadSStefano Zampini     PetscFunctionReturn(0);
13063c07aadSStefano Zampini   }
13163c07aadSStefano Zampini 
13263c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
13363c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
13463c07aadSStefano Zampini   for (i=rstart; i<rend; i++) {
13563c07aadSStefano Zampini     ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
13663c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
13763c07aadSStefano Zampini     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
13863c07aadSStefano Zampini   }
13963c07aadSStefano Zampini   PetscFunctionReturn(0);
14063c07aadSStefano Zampini }
14163c07aadSStefano Zampini 
14263c07aadSStefano Zampini #undef __FUNCT__
14363c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ"
14463c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij)
14563c07aadSStefano Zampini {
14663c07aadSStefano Zampini   PetscErrorCode        ierr;
14763c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag = (Mat_SeqAIJ*)A->data;
148*ea9daf28SStefano Zampini   int                   type;
14963c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
15063c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
15163c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag;
15263c07aadSStefano Zampini 
15363c07aadSStefano Zampini   PetscFunctionBegin;
15463c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
155*ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
156*ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
157*ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
15863c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
15963c07aadSStefano Zampini   /*
16063c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
16163c07aadSStefano Zampini   */
16263c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));
16363c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));
16463c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
165*ea9daf28SStefano Zampini 
166*ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
16763c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
16863c07aadSStefano Zampini   PetscFunctionReturn(0);
16963c07aadSStefano Zampini }
17063c07aadSStefano Zampini 
17163c07aadSStefano Zampini #undef __FUNCT__
17263c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ"
17363c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
17463c07aadSStefano Zampini {
17563c07aadSStefano Zampini   PetscErrorCode        ierr;
17663c07aadSStefano Zampini   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
17763c07aadSStefano Zampini   Mat_SeqAIJ            *pdiag,*poffd;
17863c07aadSStefano Zampini   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
179*ea9daf28SStefano Zampini   int                   type;
18063c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
18163c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
18263c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag,*hoffd;
18363c07aadSStefano Zampini 
18463c07aadSStefano Zampini   PetscFunctionBegin;
18563c07aadSStefano Zampini   pdiag = (Mat_SeqAIJ*) pA->A->data;
18663c07aadSStefano Zampini   poffd = (Mat_SeqAIJ*) pA->B->data;
18763c07aadSStefano Zampini   /* cstart is only valid for square MPIAIJ layed out in the usual way */
18863c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
18963c07aadSStefano Zampini 
19063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
191*ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
192*ea9daf28SStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
193*ea9daf28SStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
19463c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
19563c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
19663c07aadSStefano Zampini 
19763c07aadSStefano Zampini   /*
19863c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
19963c07aadSStefano Zampini   */
20063c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
20163c07aadSStefano Zampini   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
20263c07aadSStefano Zampini   jj  = (PetscInt*)hdiag->j;
20363c07aadSStefano Zampini   pjj = (PetscInt*)pdiag->j;
20463c07aadSStefano Zampini   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
20563c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
20663c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
20763c07aadSStefano Zampini   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
20863c07aadSStefano Zampini      If we hacked a hypre a bit more we might be able to avoid this step */
20963c07aadSStefano Zampini   jj  = (PetscInt*) hoffd->j;
21063c07aadSStefano Zampini   pjj = (PetscInt*) poffd->j;
21163c07aadSStefano Zampini   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
21263c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
21363c07aadSStefano Zampini 
214*ea9daf28SStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
21563c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
21663c07aadSStefano Zampini   PetscFunctionReturn(0);
21763c07aadSStefano Zampini }
21863c07aadSStefano Zampini 
21963c07aadSStefano Zampini #undef __FUNCT__
22063c07aadSStefano Zampini #define __FUNCT__ "MatConvert_AIJ_HYPRE"
22163c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
22263c07aadSStefano Zampini {
22363c07aadSStefano Zampini   Mat_HYPRE      *hB;
22463c07aadSStefano Zampini   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
22563c07aadSStefano Zampini   PetscErrorCode ierr;
22663c07aadSStefano Zampini 
22763c07aadSStefano Zampini   PetscFunctionBegin;
22863c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
22963c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
23063c07aadSStefano Zampini     /* always destroy the old matrix and create a new memory;
23163c07aadSStefano Zampini        hope this does not churn the memory too much. The problem
23263c07aadSStefano Zampini        is I do not know if it is possible to put the matrix back to
23363c07aadSStefano Zampini        its initial state so that we can directly copy the values
23463c07aadSStefano Zampini        the second time through. */
23563c07aadSStefano Zampini     hB = (Mat_HYPRE*)((*B)->data);
23663c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
23763c07aadSStefano Zampini   } else {
23863c07aadSStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
23963c07aadSStefano Zampini     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
24063c07aadSStefano Zampini     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
24163c07aadSStefano Zampini     ierr = MatSetUp(*B);CHKERRQ(ierr);
24263c07aadSStefano Zampini     hB   = (Mat_HYPRE*)((*B)->data);
24363c07aadSStefano Zampini   }
24463c07aadSStefano Zampini   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
24563c07aadSStefano Zampini   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
24663c07aadSStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
24763c07aadSStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
24863c07aadSStefano Zampini   PetscFunctionReturn(0);
24963c07aadSStefano Zampini }
25063c07aadSStefano Zampini 
25163c07aadSStefano Zampini #undef __FUNCT__
25263c07aadSStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_AIJ"
253*ea9daf28SStefano Zampini static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
25463c07aadSStefano Zampini {
25563c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
25663c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
25763c07aadSStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
25863c07aadSStefano Zampini   MPI_Comm           comm;
25963c07aadSStefano Zampini   PetscScalar        *da,*oa,*aptr;
26063c07aadSStefano Zampini   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
26163c07aadSStefano Zampini   PetscInt           i,dnnz,onnz,m,n;
26263c07aadSStefano Zampini   int                type;
26363c07aadSStefano Zampini   PetscMPIInt        size;
26463c07aadSStefano Zampini   PetscErrorCode     ierr;
26563c07aadSStefano Zampini 
26663c07aadSStefano Zampini   PetscFunctionBegin;
26763c07aadSStefano Zampini   comm = PetscObjectComm((PetscObject)A);
26863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
26963c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
27063c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
27163c07aadSStefano Zampini     PetscBool ismpiaij,isseqaij;
27263c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
27363c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
27463c07aadSStefano Zampini     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
27563c07aadSStefano Zampini   }
27663c07aadSStefano Zampini   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
27763c07aadSStefano Zampini 
27863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
27963c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(parcsr);
28063c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(parcsr);
28163c07aadSStefano Zampini   m     = hypre_CSRMatrixNumRows(hdiag);
28263c07aadSStefano Zampini   n     = hypre_CSRMatrixNumCols(hdiag);
28363c07aadSStefano Zampini   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
28463c07aadSStefano Zampini   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
28563c07aadSStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
28663c07aadSStefano Zampini     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
28763c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
28863c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
28963c07aadSStefano Zampini   } else {
29063c07aadSStefano Zampini     PetscInt  nr;
29163c07aadSStefano Zampini     PetscBool done;
29263c07aadSStefano Zampini     if (size > 1) {
29363c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
29463c07aadSStefano Zampini 
29563c07aadSStefano Zampini       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
29663c07aadSStefano 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);
29763c07aadSStefano 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);
29863c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
29963c07aadSStefano Zampini     } else {
30063c07aadSStefano Zampini       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
30163c07aadSStefano Zampini       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
30263c07aadSStefano 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);
30363c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
30463c07aadSStefano Zampini     }
30563c07aadSStefano Zampini   }
30663c07aadSStefano Zampini   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
30763c07aadSStefano Zampini   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
30863c07aadSStefano Zampini   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
30963c07aadSStefano Zampini   iptr = djj;
31063c07aadSStefano Zampini   aptr = da;
31163c07aadSStefano Zampini   for (i=0; i<m; i++) {
31263c07aadSStefano Zampini     PetscInt nc = dii[i+1]-dii[i];
31363c07aadSStefano Zampini     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
31463c07aadSStefano Zampini     iptr += nc;
31563c07aadSStefano Zampini     aptr += nc;
31663c07aadSStefano Zampini   }
31763c07aadSStefano Zampini   if (size > 1) {
31863c07aadSStefano Zampini     PetscInt *offdj,*coffd;
31963c07aadSStefano Zampini 
32063c07aadSStefano Zampini     if (reuse != MAT_REUSE_MATRIX) {
32163c07aadSStefano Zampini       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
32263c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
32363c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
32463c07aadSStefano Zampini     } else {
32563c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
32663c07aadSStefano Zampini       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
32763c07aadSStefano Zampini       PetscBool  done;
32863c07aadSStefano Zampini 
32963c07aadSStefano Zampini       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
33063c07aadSStefano 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);
33163c07aadSStefano 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);
33263c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
33363c07aadSStefano Zampini     }
33463c07aadSStefano Zampini     ierr  = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
33563c07aadSStefano Zampini     offdj = hypre_CSRMatrixJ(hoffd);
33663c07aadSStefano Zampini     coffd = hypre_ParCSRMatrixColMapOffd(parcsr);
33763c07aadSStefano Zampini     for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]];
33863c07aadSStefano Zampini     ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr);
33963c07aadSStefano Zampini     iptr = ojj;
34063c07aadSStefano Zampini     aptr = oa;
34163c07aadSStefano Zampini     for (i=0; i<m; i++) {
34263c07aadSStefano Zampini        PetscInt nc = oii[i+1]-oii[i];
34363c07aadSStefano Zampini        ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
34463c07aadSStefano Zampini        iptr += nc;
34563c07aadSStefano Zampini        aptr += nc;
34663c07aadSStefano Zampini     }
34763c07aadSStefano Zampini     if (reuse != MAT_REUSE_MATRIX) {
34863c07aadSStefano Zampini       Mat_MPIAIJ *b;
34963c07aadSStefano Zampini       Mat_SeqAIJ *d,*o;
35063c07aadSStefano Zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,
35163c07aadSStefano Zampini                                             djj,da,oii,ojj,oa,B);CHKERRQ(ierr);
35263c07aadSStefano Zampini       /* hack MPIAIJ */
35363c07aadSStefano Zampini       b          = (Mat_MPIAIJ*)((*B)->data);
35463c07aadSStefano Zampini       d          = (Mat_SeqAIJ*)b->A->data;
35563c07aadSStefano Zampini       o          = (Mat_SeqAIJ*)b->B->data;
35663c07aadSStefano Zampini       d->free_a  = PETSC_TRUE;
35763c07aadSStefano Zampini       d->free_ij = PETSC_TRUE;
35863c07aadSStefano Zampini       o->free_a  = PETSC_TRUE;
35963c07aadSStefano Zampini       o->free_ij = PETSC_TRUE;
36063c07aadSStefano Zampini     }
36163c07aadSStefano Zampini   } else if (reuse != MAT_REUSE_MATRIX) {
36263c07aadSStefano Zampini     Mat_SeqAIJ* b;
36363c07aadSStefano Zampini     ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
36463c07aadSStefano Zampini     /* hack SeqAIJ */
36563c07aadSStefano Zampini     b          = (Mat_SeqAIJ*)((*B)->data);
36663c07aadSStefano Zampini     b->free_a  = PETSC_TRUE;
36763c07aadSStefano Zampini     b->free_ij = PETSC_TRUE;
36863c07aadSStefano Zampini   }
36963c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
37063c07aadSStefano Zampini     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
37163c07aadSStefano Zampini   }
37263c07aadSStefano Zampini   PetscFunctionReturn(0);
37363c07aadSStefano Zampini }
37463c07aadSStefano Zampini 
37563c07aadSStefano Zampini #undef __FUNCT__
37663c07aadSStefano Zampini #define __FUNCT__ "MatMultTranspose_HYPRE"
377*ea9daf28SStefano Zampini static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
37863c07aadSStefano Zampini {
37963c07aadSStefano Zampini   PetscErrorCode ierr;
38063c07aadSStefano Zampini 
38163c07aadSStefano Zampini   PetscFunctionBegin;
38263c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
38363c07aadSStefano Zampini   PetscFunctionReturn(0);
38463c07aadSStefano Zampini }
38563c07aadSStefano Zampini 
38663c07aadSStefano Zampini #undef __FUNCT__
38763c07aadSStefano Zampini #define __FUNCT__ "MatMult_HYPRE"
388*ea9daf28SStefano Zampini static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
38963c07aadSStefano Zampini {
39063c07aadSStefano Zampini   PetscErrorCode ierr;
39163c07aadSStefano Zampini 
39263c07aadSStefano Zampini   PetscFunctionBegin;
39363c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
39463c07aadSStefano Zampini   PetscFunctionReturn(0);
39563c07aadSStefano Zampini }
39663c07aadSStefano Zampini 
39763c07aadSStefano Zampini #undef __FUNCT__
39863c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_MultKernel_Private"
399*ea9daf28SStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
40063c07aadSStefano Zampini {
40163c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
40263c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
40363c07aadSStefano Zampini   hypre_ParVector    *hx,*hy;
40463c07aadSStefano Zampini   PetscScalar        *ax,*ay,*sax,*say;
40563c07aadSStefano Zampini   PetscErrorCode     ierr;
40663c07aadSStefano Zampini 
40763c07aadSStefano Zampini   PetscFunctionBegin;
40863c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
40963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
41063c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
41163c07aadSStefano Zampini   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
41263c07aadSStefano Zampini   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
41363c07aadSStefano Zampini   if (trans) {
41463c07aadSStefano Zampini     HYPREReplacePointer(hA->x,ay,say);
41563c07aadSStefano Zampini     HYPREReplacePointer(hA->b,ax,sax);
41663c07aadSStefano Zampini     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
41763c07aadSStefano Zampini     HYPREReplacePointer(hA->x,say,ay);
41863c07aadSStefano Zampini     HYPREReplacePointer(hA->b,sax,ax);
41963c07aadSStefano Zampini   } else {
42063c07aadSStefano Zampini     HYPREReplacePointer(hA->x,ax,sax);
42163c07aadSStefano Zampini     HYPREReplacePointer(hA->b,ay,say);
42263c07aadSStefano Zampini     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
42363c07aadSStefano Zampini     HYPREReplacePointer(hA->x,sax,ax);
42463c07aadSStefano Zampini     HYPREReplacePointer(hA->b,say,ay);
42563c07aadSStefano Zampini   }
42663c07aadSStefano Zampini   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
42763c07aadSStefano Zampini   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
42863c07aadSStefano Zampini   PetscFunctionReturn(0);
42963c07aadSStefano Zampini }
43063c07aadSStefano Zampini 
43163c07aadSStefano Zampini #undef __FUNCT__
43263c07aadSStefano Zampini #define __FUNCT__ "MatDestroy_HYPRE"
433*ea9daf28SStefano Zampini static PetscErrorCode MatDestroy_HYPRE(Mat A)
43463c07aadSStefano Zampini {
43563c07aadSStefano Zampini   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
43663c07aadSStefano Zampini   PetscErrorCode ierr;
43763c07aadSStefano Zampini 
43863c07aadSStefano Zampini   PetscFunctionBegin;
43963c07aadSStefano Zampini   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
44063c07aadSStefano Zampini   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
44163c07aadSStefano Zampini   if (hA->ij) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
44263c07aadSStefano Zampini   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr);}
44363c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
44463c07aadSStefano Zampini   ierr = PetscFree(A->data);CHKERRQ(ierr);
44563c07aadSStefano Zampini   PetscFunctionReturn(0);
44663c07aadSStefano Zampini }
44763c07aadSStefano Zampini 
44863c07aadSStefano Zampini #undef __FUNCT__
44963c07aadSStefano Zampini #define __FUNCT__ "MatSetUp_HYPRE"
450*ea9daf28SStefano Zampini static PetscErrorCode MatSetUp_HYPRE(Mat A)
45163c07aadSStefano Zampini {
45263c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
45363c07aadSStefano Zampini   Vec                x,b;
45463c07aadSStefano Zampini   PetscErrorCode     ierr;
45563c07aadSStefano Zampini 
45663c07aadSStefano Zampini   PetscFunctionBegin;
45763c07aadSStefano Zampini   if (hA->x) PetscFunctionReturn(0);
45863c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
45963c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
46063c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
46163c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
46263c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
46363c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
46463c07aadSStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
46563c07aadSStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
46663c07aadSStefano Zampini   PetscFunctionReturn(0);
46763c07aadSStefano Zampini }
46863c07aadSStefano Zampini 
46963c07aadSStefano Zampini #undef __FUNCT__
47063c07aadSStefano Zampini #define __FUNCT__ "MatAssemblyEnd_HYPRE"
471*ea9daf28SStefano Zampini static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
47263c07aadSStefano Zampini {
47363c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
47463c07aadSStefano Zampini   PetscErrorCode     ierr;
47563c07aadSStefano Zampini 
47663c07aadSStefano Zampini   PetscFunctionBegin;
47763c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
47863c07aadSStefano Zampini   PetscFunctionReturn(0);
47963c07aadSStefano Zampini }
48063c07aadSStefano Zampini 
48163c07aadSStefano Zampini /* Is this needed?
48263c07aadSStefano Zampini   int                type;
48363c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
48463c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
48563c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR supported");
48663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
48763c07aadSStefano Zampini   PetscStackCallStandard(hypre_MatvecCommPkgCreate,(parcsr));
48863c07aadSStefano Zampini */
48963c07aadSStefano Zampini #undef __FUNCT__
49063c07aadSStefano Zampini #define __FUNCT__ "MatCreate_HYPRE"
49163c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
49263c07aadSStefano Zampini {
49363c07aadSStefano Zampini   Mat_HYPRE      *hB;
49463c07aadSStefano Zampini   PetscErrorCode ierr;
49563c07aadSStefano Zampini 
49663c07aadSStefano Zampini   PetscFunctionBegin;
49763c07aadSStefano Zampini   ierr          = PetscNewLog(B,&hB);CHKERRQ(ierr);
49863c07aadSStefano Zampini   B->data       = (void*)hB;
49963c07aadSStefano Zampini   B->rmap->bs   = 1;
50063c07aadSStefano Zampini   B->assembled  = PETSC_FALSE;
50163c07aadSStefano Zampini 
50263c07aadSStefano Zampini   B->ops->mult          = MatMult_HYPRE;
50363c07aadSStefano Zampini   B->ops->multtranspose = MatMultTranspose_HYPRE;
50463c07aadSStefano Zampini   B->ops->setup         = MatSetUp_HYPRE;
50563c07aadSStefano Zampini   B->ops->destroy       = MatDestroy_HYPRE;
50663c07aadSStefano Zampini   B->ops->assemblyend   = MatAssemblyEnd_HYPRE;
50763c07aadSStefano Zampini 
50863c07aadSStefano Zampini   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
50963c07aadSStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
51063c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
51163c07aadSStefano Zampini   PetscFunctionReturn(0);
51263c07aadSStefano Zampini }
51363c07aadSStefano Zampini 
51463c07aadSStefano Zampini #if 0
51563c07aadSStefano Zampini /*
51663c07aadSStefano Zampini     Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format
51763c07aadSStefano Zampini 
51863c07aadSStefano Zampini     This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first
51963c07aadSStefano Zampini     which will corrupt the PETSc data structure if we did this. Need a work around to this problem.
52063c07aadSStefano Zampini */
52163c07aadSStefano Zampini #include <_hypre_IJ_mv.h>
52263c07aadSStefano Zampini #include <HYPRE_IJ_mv.h>
52363c07aadSStefano Zampini 
52463c07aadSStefano Zampini #undef __FUNCT__
52563c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixLink"
52663c07aadSStefano Zampini PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij)
52763c07aadSStefano Zampini {
52863c07aadSStefano Zampini   PetscErrorCode        ierr;
52963c07aadSStefano Zampini   PetscInt              rstart,rend,cstart,cend;
53063c07aadSStefano Zampini   PetscBool             flg;
53163c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
53263c07aadSStefano Zampini 
53363c07aadSStefano Zampini   PetscFunctionBegin;
53463c07aadSStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
53563c07aadSStefano Zampini   PetscValidType(A,1);
53663c07aadSStefano Zampini   PetscValidPointer(ij,2);
53763c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
53863c07aadSStefano Zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices");
53963c07aadSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
54063c07aadSStefano Zampini 
54163c07aadSStefano Zampini   ierr   = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
54263c07aadSStefano Zampini   rstart = A->rmap->rstart;
54363c07aadSStefano Zampini   rend   = A->rmap->rend;
54463c07aadSStefano Zampini   cstart = A->cmap->rstart;
54563c07aadSStefano Zampini   cend   = A->cmap->rend;
54663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij));
54763c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
54863c07aadSStefano Zampini 
54963c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij));
55063c07aadSStefano Zampini   PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij));
55163c07aadSStefano Zampini 
55263c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
55363c07aadSStefano Zampini 
55463c07aadSStefano Zampini   /* this is the Hack part where we monkey directly with the hypre datastructures */
55563c07aadSStefano Zampini 
55663c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij));
55763c07aadSStefano Zampini   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
55863c07aadSStefano Zampini   PetscFunctionReturn(0);
55963c07aadSStefano Zampini }
56063c07aadSStefano Zampini #endif
561