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