xref: /petsc/src/mat/impls/hypre/mhypre.c (revision 63c07aad33d943fe85193412d077a1746a7c55aa)
1*63c07aadSStefano Zampini 
2*63c07aadSStefano Zampini /*
3*63c07aadSStefano Zampini     Creates hypre ijmatrix from PETSc matrix
4*63c07aadSStefano Zampini */
5*63c07aadSStefano Zampini #include <petscsys.h>
6*63c07aadSStefano Zampini #include <petsc/private/matimpl.h>
7*63c07aadSStefano Zampini #include <../src/mat/impls/hypre/mhypre.h>
8*63c07aadSStefano Zampini #include <../src/mat/impls/aij/mpi/mpiaij.h>
9*63c07aadSStefano Zampini #include <HYPRE_struct_mv.h>
10*63c07aadSStefano Zampini #include <HYPRE_struct_ls.h>
11*63c07aadSStefano Zampini #include <_hypre_struct_mv.h>
12*63c07aadSStefano Zampini 
13*63c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*);
14*63c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix);
15*63c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
16*63c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
17*63c07aadSStefano Zampini static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool);
18*63c07aadSStefano Zampini 
19*63c07aadSStefano Zampini #undef __FUNCT__
20*63c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate"
21*63c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij)
22*63c07aadSStefano Zampini {
23*63c07aadSStefano Zampini   PetscErrorCode ierr;
24*63c07aadSStefano Zampini   PetscInt       i,n_d,n_o;
25*63c07aadSStefano Zampini   const PetscInt *ia_d,*ia_o;
26*63c07aadSStefano Zampini   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
27*63c07aadSStefano Zampini   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
28*63c07aadSStefano Zampini 
29*63c07aadSStefano Zampini   PetscFunctionBegin;
30*63c07aadSStefano Zampini   if (A_d) { /* determine number of nonzero entries in local diagonal part */
31*63c07aadSStefano Zampini     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
32*63c07aadSStefano Zampini     if (done_d) {
33*63c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr);
34*63c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
35*63c07aadSStefano Zampini         nnz_d[i] = ia_d[i+1] - ia_d[i];
36*63c07aadSStefano Zampini       }
37*63c07aadSStefano Zampini     }
38*63c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
39*63c07aadSStefano Zampini   }
40*63c07aadSStefano Zampini   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
41*63c07aadSStefano Zampini     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
42*63c07aadSStefano Zampini     if (done_o) {
43*63c07aadSStefano Zampini       ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr);
44*63c07aadSStefano Zampini       for (i=0; i<n_o; i++) {
45*63c07aadSStefano Zampini         nnz_o[i] = ia_o[i+1] - ia_o[i];
46*63c07aadSStefano Zampini       }
47*63c07aadSStefano Zampini     }
48*63c07aadSStefano Zampini     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
49*63c07aadSStefano Zampini   }
50*63c07aadSStefano Zampini   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
51*63c07aadSStefano Zampini     if (!done_o) { /* only diagonal part */
52*63c07aadSStefano Zampini       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
53*63c07aadSStefano Zampini       for (i=0; i<n_d; i++) {
54*63c07aadSStefano Zampini         nnz_o[i] = 0;
55*63c07aadSStefano Zampini       }
56*63c07aadSStefano Zampini     }
57*63c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o));
58*63c07aadSStefano Zampini     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
59*63c07aadSStefano Zampini     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
60*63c07aadSStefano Zampini   }
61*63c07aadSStefano Zampini   PetscFunctionReturn(0);
62*63c07aadSStefano Zampini }
63*63c07aadSStefano Zampini 
64*63c07aadSStefano Zampini #undef __FUNCT__
65*63c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_CreateFromMat"
66*63c07aadSStefano Zampini static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA)
67*63c07aadSStefano Zampini {
68*63c07aadSStefano Zampini   PetscErrorCode ierr;
69*63c07aadSStefano Zampini   PetscInt       rstart,rend,cstart,cend;
70*63c07aadSStefano Zampini 
71*63c07aadSStefano Zampini   PetscFunctionBegin;
72*63c07aadSStefano Zampini   ierr   = MatSetUp(A);CHKERRQ(ierr);
73*63c07aadSStefano Zampini   rstart = A->rmap->rstart;
74*63c07aadSStefano Zampini   rend   = A->rmap->rend;
75*63c07aadSStefano Zampini   cstart = A->cmap->rstart;
76*63c07aadSStefano Zampini   cend   = A->cmap->rend;
77*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
78*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
79*63c07aadSStefano Zampini   {
80*63c07aadSStefano Zampini     PetscBool      same;
81*63c07aadSStefano Zampini     Mat            A_d,A_o;
82*63c07aadSStefano Zampini     const PetscInt *colmap;
83*63c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
84*63c07aadSStefano Zampini     if (same) {
85*63c07aadSStefano Zampini       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
86*63c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
87*63c07aadSStefano Zampini       PetscFunctionReturn(0);
88*63c07aadSStefano Zampini     }
89*63c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
90*63c07aadSStefano Zampini     if (same) {
91*63c07aadSStefano Zampini       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
92*63c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
93*63c07aadSStefano Zampini       PetscFunctionReturn(0);
94*63c07aadSStefano Zampini     }
95*63c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
96*63c07aadSStefano Zampini     if (same) {
97*63c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
98*63c07aadSStefano Zampini       PetscFunctionReturn(0);
99*63c07aadSStefano Zampini     }
100*63c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
101*63c07aadSStefano Zampini     if (same) {
102*63c07aadSStefano Zampini       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
103*63c07aadSStefano Zampini       PetscFunctionReturn(0);
104*63c07aadSStefano Zampini     }
105*63c07aadSStefano Zampini   }
106*63c07aadSStefano Zampini   PetscFunctionReturn(0);
107*63c07aadSStefano Zampini }
108*63c07aadSStefano Zampini 
109*63c07aadSStefano Zampini #undef __FUNCT__
110*63c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixCopy"
111*63c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij)
112*63c07aadSStefano Zampini {
113*63c07aadSStefano Zampini   PetscErrorCode    ierr;
114*63c07aadSStefano Zampini   PetscInt          i,rstart,rend,ncols,nr,nc;
115*63c07aadSStefano Zampini   const PetscScalar *values;
116*63c07aadSStefano Zampini   const PetscInt    *cols;
117*63c07aadSStefano Zampini   PetscBool         flg;
118*63c07aadSStefano Zampini 
119*63c07aadSStefano Zampini   PetscFunctionBegin;
120*63c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
121*63c07aadSStefano Zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
122*63c07aadSStefano Zampini   if (flg && nr == nc) {
123*63c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
124*63c07aadSStefano Zampini     PetscFunctionReturn(0);
125*63c07aadSStefano Zampini   }
126*63c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
127*63c07aadSStefano Zampini   if (flg) {
128*63c07aadSStefano Zampini     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
129*63c07aadSStefano Zampini     PetscFunctionReturn(0);
130*63c07aadSStefano Zampini   }
131*63c07aadSStefano Zampini 
132*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
133*63c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
134*63c07aadSStefano Zampini   for (i=rstart; i<rend; i++) {
135*63c07aadSStefano Zampini     ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
136*63c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
137*63c07aadSStefano Zampini     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
138*63c07aadSStefano Zampini   }
139*63c07aadSStefano Zampini   PetscFunctionReturn(0);
140*63c07aadSStefano Zampini }
141*63c07aadSStefano Zampini 
142*63c07aadSStefano Zampini #undef __FUNCT__
143*63c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ"
144*63c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij)
145*63c07aadSStefano Zampini {
146*63c07aadSStefano Zampini   PetscErrorCode ierr;
147*63c07aadSStefano Zampini   Mat_SeqAIJ     *pdiag = (Mat_SeqAIJ*)A->data;
148*63c07aadSStefano Zampini 
149*63c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
150*63c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
151*63c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag;
152*63c07aadSStefano Zampini 
153*63c07aadSStefano Zampini   PetscFunctionBegin;
154*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
155*63c07aadSStefano Zampini   par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij);
156*63c07aadSStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
157*63c07aadSStefano Zampini   hdiag      = hypre_ParCSRMatrixDiag(par_matrix);
158*63c07aadSStefano Zampini   /*
159*63c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
160*63c07aadSStefano Zampini   */
161*63c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));
162*63c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));
163*63c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
164*63c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
165*63c07aadSStefano Zampini   PetscFunctionReturn(0);
166*63c07aadSStefano Zampini }
167*63c07aadSStefano Zampini 
168*63c07aadSStefano Zampini #undef __FUNCT__
169*63c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ"
170*63c07aadSStefano Zampini static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
171*63c07aadSStefano Zampini {
172*63c07aadSStefano Zampini   PetscErrorCode ierr;
173*63c07aadSStefano Zampini   Mat_MPIAIJ     *pA = (Mat_MPIAIJ*)A->data;
174*63c07aadSStefano Zampini   Mat_SeqAIJ     *pdiag,*poffd;
175*63c07aadSStefano Zampini   PetscInt       i,*garray = pA->garray,*jj,cstart,*pjj;
176*63c07aadSStefano Zampini 
177*63c07aadSStefano Zampini   hypre_ParCSRMatrix    *par_matrix;
178*63c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
179*63c07aadSStefano Zampini   hypre_CSRMatrix       *hdiag,*hoffd;
180*63c07aadSStefano Zampini 
181*63c07aadSStefano Zampini   PetscFunctionBegin;
182*63c07aadSStefano Zampini   pdiag = (Mat_SeqAIJ*) pA->A->data;
183*63c07aadSStefano Zampini   poffd = (Mat_SeqAIJ*) pA->B->data;
184*63c07aadSStefano Zampini   /* cstart is only valid for square MPIAIJ layed out in the usual way */
185*63c07aadSStefano Zampini   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
186*63c07aadSStefano Zampini 
187*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
188*63c07aadSStefano Zampini   par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij);
189*63c07aadSStefano Zampini   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
190*63c07aadSStefano Zampini   hdiag      = hypre_ParCSRMatrixDiag(par_matrix);
191*63c07aadSStefano Zampini   hoffd      = hypre_ParCSRMatrixOffd(par_matrix);
192*63c07aadSStefano Zampini 
193*63c07aadSStefano Zampini   /*
194*63c07aadSStefano Zampini        this is the Hack part where we monkey directly with the hypre datastructures
195*63c07aadSStefano Zampini   */
196*63c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
197*63c07aadSStefano Zampini   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
198*63c07aadSStefano Zampini   jj  = (PetscInt*)hdiag->j;
199*63c07aadSStefano Zampini   pjj = (PetscInt*)pdiag->j;
200*63c07aadSStefano Zampini   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
201*63c07aadSStefano Zampini   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
202*63c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
203*63c07aadSStefano Zampini   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
204*63c07aadSStefano Zampini      If we hacked a hypre a bit more we might be able to avoid this step */
205*63c07aadSStefano Zampini   jj  = (PetscInt*) hoffd->j;
206*63c07aadSStefano Zampini   pjj = (PetscInt*) poffd->j;
207*63c07aadSStefano Zampini   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
208*63c07aadSStefano Zampini   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
209*63c07aadSStefano Zampini 
210*63c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
211*63c07aadSStefano Zampini   PetscFunctionReturn(0);
212*63c07aadSStefano Zampini }
213*63c07aadSStefano Zampini 
214*63c07aadSStefano Zampini #undef __FUNCT__
215*63c07aadSStefano Zampini #define __FUNCT__ "MatConvert_AIJ_HYPRE"
216*63c07aadSStefano Zampini PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
217*63c07aadSStefano Zampini {
218*63c07aadSStefano Zampini   Mat_HYPRE      *hB;
219*63c07aadSStefano Zampini   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
220*63c07aadSStefano Zampini   PetscErrorCode ierr;
221*63c07aadSStefano Zampini 
222*63c07aadSStefano Zampini   PetscFunctionBegin;
223*63c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
224*63c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
225*63c07aadSStefano Zampini     /* always destroy the old matrix and create a new memory;
226*63c07aadSStefano Zampini        hope this does not churn the memory too much. The problem
227*63c07aadSStefano Zampini        is I do not know if it is possible to put the matrix back to
228*63c07aadSStefano Zampini        its initial state so that we can directly copy the values
229*63c07aadSStefano Zampini        the second time through. */
230*63c07aadSStefano Zampini     hB = (Mat_HYPRE*)((*B)->data);
231*63c07aadSStefano Zampini     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
232*63c07aadSStefano Zampini   } else {
233*63c07aadSStefano Zampini     ierr = MatCreate(comm,B);CHKERRQ(ierr);
234*63c07aadSStefano Zampini     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
235*63c07aadSStefano Zampini     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
236*63c07aadSStefano Zampini     ierr = MatSetUp(*B);CHKERRQ(ierr);
237*63c07aadSStefano Zampini     hB   = (Mat_HYPRE*)((*B)->data);
238*63c07aadSStefano Zampini   }
239*63c07aadSStefano Zampini   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
240*63c07aadSStefano Zampini   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
241*63c07aadSStefano Zampini   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
242*63c07aadSStefano Zampini   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
243*63c07aadSStefano Zampini   PetscFunctionReturn(0);
244*63c07aadSStefano Zampini }
245*63c07aadSStefano Zampini 
246*63c07aadSStefano Zampini #undef __FUNCT__
247*63c07aadSStefano Zampini #define __FUNCT__ "MatConvert_HYPRE_AIJ"
248*63c07aadSStefano Zampini PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
249*63c07aadSStefano Zampini {
250*63c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
251*63c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
252*63c07aadSStefano Zampini   hypre_CSRMatrix    *hdiag,*hoffd;
253*63c07aadSStefano Zampini   MPI_Comm           comm;
254*63c07aadSStefano Zampini   PetscScalar        *da,*oa,*aptr;
255*63c07aadSStefano Zampini   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
256*63c07aadSStefano Zampini   PetscInt           i,dnnz,onnz,m,n;
257*63c07aadSStefano Zampini   int                type;
258*63c07aadSStefano Zampini   PetscMPIInt        size;
259*63c07aadSStefano Zampini   PetscErrorCode     ierr;
260*63c07aadSStefano Zampini 
261*63c07aadSStefano Zampini   PetscFunctionBegin;
262*63c07aadSStefano Zampini   comm = PetscObjectComm((PetscObject)A);
263*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
264*63c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
265*63c07aadSStefano Zampini   if (reuse == MAT_REUSE_MATRIX) {
266*63c07aadSStefano Zampini     PetscBool ismpiaij,isseqaij;
267*63c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
268*63c07aadSStefano Zampini     ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
269*63c07aadSStefano Zampini     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
270*63c07aadSStefano Zampini   }
271*63c07aadSStefano Zampini   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
272*63c07aadSStefano Zampini 
273*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
274*63c07aadSStefano Zampini   hdiag = hypre_ParCSRMatrixDiag(parcsr);
275*63c07aadSStefano Zampini   hoffd = hypre_ParCSRMatrixOffd(parcsr);
276*63c07aadSStefano Zampini   m     = hypre_CSRMatrixNumRows(hdiag);
277*63c07aadSStefano Zampini   n     = hypre_CSRMatrixNumCols(hdiag);
278*63c07aadSStefano Zampini   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
279*63c07aadSStefano Zampini   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
280*63c07aadSStefano Zampini   if (reuse != MAT_REUSE_MATRIX) {
281*63c07aadSStefano Zampini     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
282*63c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
283*63c07aadSStefano Zampini     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
284*63c07aadSStefano Zampini   } else {
285*63c07aadSStefano Zampini     PetscInt  nr;
286*63c07aadSStefano Zampini     PetscBool done;
287*63c07aadSStefano Zampini     if (size > 1) {
288*63c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
289*63c07aadSStefano Zampini 
290*63c07aadSStefano Zampini       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
291*63c07aadSStefano 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);
292*63c07aadSStefano 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);
293*63c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
294*63c07aadSStefano Zampini     } else {
295*63c07aadSStefano Zampini       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
296*63c07aadSStefano Zampini       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
297*63c07aadSStefano 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);
298*63c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
299*63c07aadSStefano Zampini     }
300*63c07aadSStefano Zampini   }
301*63c07aadSStefano Zampini   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
302*63c07aadSStefano Zampini   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
303*63c07aadSStefano Zampini   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
304*63c07aadSStefano Zampini   iptr = djj;
305*63c07aadSStefano Zampini   aptr = da;
306*63c07aadSStefano Zampini   for (i=0; i<m; i++) {
307*63c07aadSStefano Zampini     PetscInt nc = dii[i+1]-dii[i];
308*63c07aadSStefano Zampini     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
309*63c07aadSStefano Zampini     iptr += nc;
310*63c07aadSStefano Zampini     aptr += nc;
311*63c07aadSStefano Zampini   }
312*63c07aadSStefano Zampini   if (size > 1) {
313*63c07aadSStefano Zampini     PetscInt *offdj,*coffd;
314*63c07aadSStefano Zampini 
315*63c07aadSStefano Zampini     if (reuse != MAT_REUSE_MATRIX) {
316*63c07aadSStefano Zampini       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
317*63c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
318*63c07aadSStefano Zampini       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
319*63c07aadSStefano Zampini     } else {
320*63c07aadSStefano Zampini       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
321*63c07aadSStefano Zampini       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
322*63c07aadSStefano Zampini       PetscBool  done;
323*63c07aadSStefano Zampini 
324*63c07aadSStefano Zampini       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
325*63c07aadSStefano 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);
326*63c07aadSStefano 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);
327*63c07aadSStefano Zampini       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
328*63c07aadSStefano Zampini     }
329*63c07aadSStefano Zampini     ierr  = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
330*63c07aadSStefano Zampini     offdj = hypre_CSRMatrixJ(hoffd);
331*63c07aadSStefano Zampini     coffd = hypre_ParCSRMatrixColMapOffd(parcsr);
332*63c07aadSStefano Zampini     for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]];
333*63c07aadSStefano Zampini     ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr);
334*63c07aadSStefano Zampini     iptr = ojj;
335*63c07aadSStefano Zampini     aptr = oa;
336*63c07aadSStefano Zampini     for (i=0; i<m; i++) {
337*63c07aadSStefano Zampini        PetscInt nc = oii[i+1]-oii[i];
338*63c07aadSStefano Zampini        ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
339*63c07aadSStefano Zampini        iptr += nc;
340*63c07aadSStefano Zampini        aptr += nc;
341*63c07aadSStefano Zampini     }
342*63c07aadSStefano Zampini     if (reuse != MAT_REUSE_MATRIX) {
343*63c07aadSStefano Zampini       Mat_MPIAIJ *b;
344*63c07aadSStefano Zampini       Mat_SeqAIJ *d,*o;
345*63c07aadSStefano Zampini       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,
346*63c07aadSStefano Zampini                                             djj,da,oii,ojj,oa,B);CHKERRQ(ierr);
347*63c07aadSStefano Zampini       /* hack MPIAIJ */
348*63c07aadSStefano Zampini       b          = (Mat_MPIAIJ*)((*B)->data);
349*63c07aadSStefano Zampini       d          = (Mat_SeqAIJ*)b->A->data;
350*63c07aadSStefano Zampini       o          = (Mat_SeqAIJ*)b->B->data;
351*63c07aadSStefano Zampini       d->free_a  = PETSC_TRUE;
352*63c07aadSStefano Zampini       d->free_ij = PETSC_TRUE;
353*63c07aadSStefano Zampini       o->free_a  = PETSC_TRUE;
354*63c07aadSStefano Zampini       o->free_ij = PETSC_TRUE;
355*63c07aadSStefano Zampini     }
356*63c07aadSStefano Zampini   } else if (reuse != MAT_REUSE_MATRIX) {
357*63c07aadSStefano Zampini     Mat_SeqAIJ* b;
358*63c07aadSStefano Zampini     ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
359*63c07aadSStefano Zampini     /* hack SeqAIJ */
360*63c07aadSStefano Zampini     b          = (Mat_SeqAIJ*)((*B)->data);
361*63c07aadSStefano Zampini     b->free_a  = PETSC_TRUE;
362*63c07aadSStefano Zampini     b->free_ij = PETSC_TRUE;
363*63c07aadSStefano Zampini   }
364*63c07aadSStefano Zampini   if (reuse == MAT_INPLACE_MATRIX) {
365*63c07aadSStefano Zampini     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
366*63c07aadSStefano Zampini   }
367*63c07aadSStefano Zampini   PetscFunctionReturn(0);
368*63c07aadSStefano Zampini }
369*63c07aadSStefano Zampini 
370*63c07aadSStefano Zampini #undef __FUNCT__
371*63c07aadSStefano Zampini #define __FUNCT__ "MatMultTranspose_HYPRE"
372*63c07aadSStefano Zampini PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
373*63c07aadSStefano Zampini {
374*63c07aadSStefano Zampini   PetscErrorCode ierr;
375*63c07aadSStefano Zampini 
376*63c07aadSStefano Zampini   PetscFunctionBegin;
377*63c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
378*63c07aadSStefano Zampini   PetscFunctionReturn(0);
379*63c07aadSStefano Zampini }
380*63c07aadSStefano Zampini 
381*63c07aadSStefano Zampini #undef __FUNCT__
382*63c07aadSStefano Zampini #define __FUNCT__ "MatMult_HYPRE"
383*63c07aadSStefano Zampini PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
384*63c07aadSStefano Zampini {
385*63c07aadSStefano Zampini   PetscErrorCode ierr;
386*63c07aadSStefano Zampini 
387*63c07aadSStefano Zampini   PetscFunctionBegin;
388*63c07aadSStefano Zampini   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
389*63c07aadSStefano Zampini   PetscFunctionReturn(0);
390*63c07aadSStefano Zampini }
391*63c07aadSStefano Zampini 
392*63c07aadSStefano Zampini #undef __FUNCT__
393*63c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_MultKernel_Private"
394*63c07aadSStefano Zampini PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
395*63c07aadSStefano Zampini {
396*63c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
397*63c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
398*63c07aadSStefano Zampini   hypre_ParVector    *hx,*hy;
399*63c07aadSStefano Zampini   PetscScalar        *ax,*ay,*sax,*say;
400*63c07aadSStefano Zampini   PetscErrorCode     ierr;
401*63c07aadSStefano Zampini 
402*63c07aadSStefano Zampini   PetscFunctionBegin;
403*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
404*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
405*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
406*63c07aadSStefano Zampini   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
407*63c07aadSStefano Zampini   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
408*63c07aadSStefano Zampini   if (trans) {
409*63c07aadSStefano Zampini     HYPREReplacePointer(hA->x,ay,say);
410*63c07aadSStefano Zampini     HYPREReplacePointer(hA->b,ax,sax);
411*63c07aadSStefano Zampini     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
412*63c07aadSStefano Zampini     HYPREReplacePointer(hA->x,say,ay);
413*63c07aadSStefano Zampini     HYPREReplacePointer(hA->b,sax,ax);
414*63c07aadSStefano Zampini   } else {
415*63c07aadSStefano Zampini     HYPREReplacePointer(hA->x,ax,sax);
416*63c07aadSStefano Zampini     HYPREReplacePointer(hA->b,ay,say);
417*63c07aadSStefano Zampini     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
418*63c07aadSStefano Zampini     HYPREReplacePointer(hA->x,sax,ax);
419*63c07aadSStefano Zampini     HYPREReplacePointer(hA->b,say,ay);
420*63c07aadSStefano Zampini   }
421*63c07aadSStefano Zampini   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
422*63c07aadSStefano Zampini   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
423*63c07aadSStefano Zampini   PetscFunctionReturn(0);
424*63c07aadSStefano Zampini }
425*63c07aadSStefano Zampini 
426*63c07aadSStefano Zampini #undef __FUNCT__
427*63c07aadSStefano Zampini #define __FUNCT__ "MatDestroy_HYPRE"
428*63c07aadSStefano Zampini PetscErrorCode MatDestroy_HYPRE(Mat A)
429*63c07aadSStefano Zampini {
430*63c07aadSStefano Zampini   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
431*63c07aadSStefano Zampini   PetscErrorCode ierr;
432*63c07aadSStefano Zampini 
433*63c07aadSStefano Zampini   PetscFunctionBegin;
434*63c07aadSStefano Zampini   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
435*63c07aadSStefano Zampini   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
436*63c07aadSStefano Zampini   if (hA->ij) PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
437*63c07aadSStefano Zampini   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr);}
438*63c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
439*63c07aadSStefano Zampini   ierr = PetscFree(A->data);CHKERRQ(ierr);
440*63c07aadSStefano Zampini   PetscFunctionReturn(0);
441*63c07aadSStefano Zampini }
442*63c07aadSStefano Zampini 
443*63c07aadSStefano Zampini #undef __FUNCT__
444*63c07aadSStefano Zampini #define __FUNCT__ "MatSetUp_HYPRE"
445*63c07aadSStefano Zampini PetscErrorCode MatSetUp_HYPRE(Mat A)
446*63c07aadSStefano Zampini {
447*63c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
448*63c07aadSStefano Zampini   Vec                x,b;
449*63c07aadSStefano Zampini   PetscErrorCode     ierr;
450*63c07aadSStefano Zampini 
451*63c07aadSStefano Zampini   PetscFunctionBegin;
452*63c07aadSStefano Zampini   if (hA->x) PetscFunctionReturn(0);
453*63c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
454*63c07aadSStefano Zampini   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
455*63c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
456*63c07aadSStefano Zampini   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
457*63c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
458*63c07aadSStefano Zampini   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
459*63c07aadSStefano Zampini   ierr = VecDestroy(&x);CHKERRQ(ierr);
460*63c07aadSStefano Zampini   ierr = VecDestroy(&b);CHKERRQ(ierr);
461*63c07aadSStefano Zampini   PetscFunctionReturn(0);
462*63c07aadSStefano Zampini }
463*63c07aadSStefano Zampini 
464*63c07aadSStefano Zampini #undef __FUNCT__
465*63c07aadSStefano Zampini #define __FUNCT__ "MatAssemblyEnd_HYPRE"
466*63c07aadSStefano Zampini PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
467*63c07aadSStefano Zampini {
468*63c07aadSStefano Zampini   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
469*63c07aadSStefano Zampini   PetscErrorCode     ierr;
470*63c07aadSStefano Zampini 
471*63c07aadSStefano Zampini   PetscFunctionBegin;
472*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
473*63c07aadSStefano Zampini   PetscFunctionReturn(0);
474*63c07aadSStefano Zampini }
475*63c07aadSStefano Zampini 
476*63c07aadSStefano Zampini /* Is this needed?
477*63c07aadSStefano Zampini   int                type;
478*63c07aadSStefano Zampini   hypre_ParCSRMatrix *parcsr;
479*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
480*63c07aadSStefano Zampini   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR supported");
481*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
482*63c07aadSStefano Zampini   PetscStackCallStandard(hypre_MatvecCommPkgCreate,(parcsr));
483*63c07aadSStefano Zampini */
484*63c07aadSStefano Zampini #undef __FUNCT__
485*63c07aadSStefano Zampini #define __FUNCT__ "MatCreate_HYPRE"
486*63c07aadSStefano Zampini PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
487*63c07aadSStefano Zampini {
488*63c07aadSStefano Zampini   Mat_HYPRE      *hB;
489*63c07aadSStefano Zampini   PetscErrorCode ierr;
490*63c07aadSStefano Zampini 
491*63c07aadSStefano Zampini   PetscFunctionBegin;
492*63c07aadSStefano Zampini   ierr          = PetscNewLog(B,&hB);CHKERRQ(ierr);
493*63c07aadSStefano Zampini   B->data       = (void*)hB;
494*63c07aadSStefano Zampini   B->rmap->bs   = 1;
495*63c07aadSStefano Zampini   B->assembled  = PETSC_FALSE;
496*63c07aadSStefano Zampini 
497*63c07aadSStefano Zampini   B->ops->mult          = MatMult_HYPRE;
498*63c07aadSStefano Zampini   B->ops->multtranspose = MatMultTranspose_HYPRE;
499*63c07aadSStefano Zampini   B->ops->setup         = MatSetUp_HYPRE;
500*63c07aadSStefano Zampini   B->ops->destroy       = MatDestroy_HYPRE;
501*63c07aadSStefano Zampini   B->ops->assemblyend   = MatAssemblyEnd_HYPRE;
502*63c07aadSStefano Zampini 
503*63c07aadSStefano Zampini   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
504*63c07aadSStefano Zampini   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
505*63c07aadSStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
506*63c07aadSStefano Zampini   PetscFunctionReturn(0);
507*63c07aadSStefano Zampini }
508*63c07aadSStefano Zampini 
509*63c07aadSStefano Zampini #if 0
510*63c07aadSStefano Zampini /*
511*63c07aadSStefano Zampini     Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format
512*63c07aadSStefano Zampini 
513*63c07aadSStefano Zampini     This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first
514*63c07aadSStefano Zampini     which will corrupt the PETSc data structure if we did this. Need a work around to this problem.
515*63c07aadSStefano Zampini */
516*63c07aadSStefano Zampini #include <_hypre_IJ_mv.h>
517*63c07aadSStefano Zampini #include <HYPRE_IJ_mv.h>
518*63c07aadSStefano Zampini 
519*63c07aadSStefano Zampini #undef __FUNCT__
520*63c07aadSStefano Zampini #define __FUNCT__ "MatHYPRE_IJMatrixLink"
521*63c07aadSStefano Zampini PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij)
522*63c07aadSStefano Zampini {
523*63c07aadSStefano Zampini   PetscErrorCode        ierr;
524*63c07aadSStefano Zampini   PetscInt              rstart,rend,cstart,cend;
525*63c07aadSStefano Zampini   PetscBool             flg;
526*63c07aadSStefano Zampini   hypre_AuxParCSRMatrix *aux_matrix;
527*63c07aadSStefano Zampini 
528*63c07aadSStefano Zampini   PetscFunctionBegin;
529*63c07aadSStefano Zampini   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
530*63c07aadSStefano Zampini   PetscValidType(A,1);
531*63c07aadSStefano Zampini   PetscValidPointer(ij,2);
532*63c07aadSStefano Zampini   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
533*63c07aadSStefano Zampini   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices");
534*63c07aadSStefano Zampini   ierr = MatSetUp(A);CHKERRQ(ierr);
535*63c07aadSStefano Zampini 
536*63c07aadSStefano Zampini   ierr   = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
537*63c07aadSStefano Zampini   rstart = A->rmap->rstart;
538*63c07aadSStefano Zampini   rend   = A->rmap->rend;
539*63c07aadSStefano Zampini   cstart = A->cmap->rstart;
540*63c07aadSStefano Zampini   cend   = A->cmap->rend;
541*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij));
542*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
543*63c07aadSStefano Zampini 
544*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij));
545*63c07aadSStefano Zampini   PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij));
546*63c07aadSStefano Zampini 
547*63c07aadSStefano Zampini   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
548*63c07aadSStefano Zampini 
549*63c07aadSStefano Zampini   /* this is the Hack part where we monkey directly with the hypre datastructures */
550*63c07aadSStefano Zampini 
551*63c07aadSStefano Zampini   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij));
552*63c07aadSStefano Zampini   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
553*63c07aadSStefano Zampini   PetscFunctionReturn(0);
554*63c07aadSStefano Zampini }
555*63c07aadSStefano Zampini #endif
556