xref: /petsc/src/mat/impls/hypre/mhypre.c (revision 613e5ff002909a1ed211d36f0183738222529a50)
1 
2 /*
3     Creates hypre ijmatrix from PETSc matrix
4 */
5 
6 /*MC
7    MATHYPRE - MATHYPRE = "hypre" - A matrix type to be used for sequential and parallel sparse matrices
8           based on the hypre IJ interface.
9 
10    Level: intermediate
11 
12 .seealso: MatCreate()
13 M*/
14 
15 #include <petscmathypre.h>
16 #include <petsc/private/matimpl.h>
17 #include <../src/mat/impls/hypre/mhypre.h>
18 #include <../src/mat/impls/aij/mpi/mpiaij.h>
19 #include <../src/vec/vec/impls/hypre/vhyp.h>
20 #include <HYPRE.h>
21 #include <HYPRE_utilities.h>
22 #include <_hypre_parcsr_ls.h>
23 
24 static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*);
25 static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix);
26 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
27 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
28 static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool);
29 static PetscErrorCode hypre_array_destroy(void*);
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate"
33 static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij)
34 {
35   PetscErrorCode ierr;
36   PetscInt       i,n_d,n_o;
37   const PetscInt *ia_d,*ia_o;
38   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
39   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
40 
41   PetscFunctionBegin;
42   if (A_d) { /* determine number of nonzero entries in local diagonal part */
43     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
44     if (done_d) {
45       ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr);
46       for (i=0; i<n_d; i++) {
47         nnz_d[i] = ia_d[i+1] - ia_d[i];
48       }
49     }
50     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
51   }
52   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
53     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
54     if (done_o) {
55       ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr);
56       for (i=0; i<n_o; i++) {
57         nnz_o[i] = ia_o[i+1] - ia_o[i];
58       }
59     }
60     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
61   }
62   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
63     if (!done_o) { /* only diagonal part */
64       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
65       for (i=0; i<n_d; i++) {
66         nnz_o[i] = 0;
67       }
68     }
69     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o));
70     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
71     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
72   }
73   PetscFunctionReturn(0);
74 }
75 
76 #undef __FUNCT__
77 #define __FUNCT__ "MatHYPRE_CreateFromMat"
78 static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA)
79 {
80   PetscErrorCode ierr;
81   PetscInt       rstart,rend,cstart,cend;
82 
83   PetscFunctionBegin;
84   ierr   = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
85   ierr   = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
86   rstart = A->rmap->rstart;
87   rend   = A->rmap->rend;
88   cstart = A->cmap->rstart;
89   cend   = A->cmap->rend;
90   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
91   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
92   {
93     PetscBool      same;
94     Mat            A_d,A_o;
95     const PetscInt *colmap;
96     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
97     if (same) {
98       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
99       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
100       PetscFunctionReturn(0);
101     }
102     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
103     if (same) {
104       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
105       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
106       PetscFunctionReturn(0);
107     }
108     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
109     if (same) {
110       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
111       PetscFunctionReturn(0);
112     }
113     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
114     if (same) {
115       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
116       PetscFunctionReturn(0);
117     }
118   }
119   PetscFunctionReturn(0);
120 }
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "MatHYPRE_IJMatrixCopy"
124 static PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A, HYPRE_IJMatrix ij)
125 {
126   PetscErrorCode    ierr;
127   PetscInt          i,rstart,rend,ncols,nr,nc;
128   const PetscScalar *values;
129   const PetscInt    *cols;
130   PetscBool         flg;
131 
132   PetscFunctionBegin;
133   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
134   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
135   if (flg && nr == nc) {
136     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
137     PetscFunctionReturn(0);
138   }
139   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
140   if (flg) {
141     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
142     PetscFunctionReturn(0);
143   }
144 
145   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
146   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
147   for (i=rstart; i<rend; i++) {
148     ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
149     PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
150     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
151   }
152   PetscFunctionReturn(0);
153 }
154 
155 #undef __FUNCT__
156 #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ"
157 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij)
158 {
159   PetscErrorCode        ierr;
160   Mat_SeqAIJ            *pdiag = (Mat_SeqAIJ*)A->data;
161   HYPRE_Int             type;
162   hypre_ParCSRMatrix    *par_matrix;
163   hypre_AuxParCSRMatrix *aux_matrix;
164   hypre_CSRMatrix       *hdiag;
165 
166   PetscFunctionBegin;
167   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
168   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
169   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
170   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
171   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
172   /*
173        this is the Hack part where we monkey directly with the hypre datastructures
174   */
175   ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));
176   ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));
177   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
178 
179   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
180   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
181   PetscFunctionReturn(0);
182 }
183 
184 #undef __FUNCT__
185 #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ"
186 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
187 {
188   PetscErrorCode        ierr;
189   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
190   Mat_SeqAIJ            *pdiag,*poffd;
191   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
192   HYPRE_Int             type;
193   hypre_ParCSRMatrix    *par_matrix;
194   hypre_AuxParCSRMatrix *aux_matrix;
195   hypre_CSRMatrix       *hdiag,*hoffd;
196 
197   PetscFunctionBegin;
198   pdiag = (Mat_SeqAIJ*) pA->A->data;
199   poffd = (Mat_SeqAIJ*) pA->B->data;
200   /* cstart is only valid for square MPIAIJ layed out in the usual way */
201   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
202 
203   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
204   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
205   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
206   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
207   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
208   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
209 
210   /*
211        this is the Hack part where we monkey directly with the hypre datastructures
212   */
213   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
214   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
215   jj  = (PetscInt*)hdiag->j;
216   pjj = (PetscInt*)pdiag->j;
217   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
218   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
219   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
220   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
221      If we hacked a hypre a bit more we might be able to avoid this step */
222   jj  = (PetscInt*) hoffd->j;
223   pjj = (PetscInt*) poffd->j;
224   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
225   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
226 
227   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
228   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
229   PetscFunctionReturn(0);
230 }
231 
232 #undef __FUNCT__
233 #define __FUNCT__ "MatConvert_HYPRE_IS"
234 static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B)
235 {
236   Mat_HYPRE*             mhA = (Mat_HYPRE*)(A->data);
237   Mat                    lA;
238   ISLocalToGlobalMapping rl2g,cl2g;
239   IS                     is;
240   hypre_ParCSRMatrix     *hA;
241   hypre_CSRMatrix        *hdiag,*hoffd;
242   MPI_Comm               comm;
243   PetscScalar            *hdd,*hod,*aa,*data;
244   HYPRE_Int              *col_map_offd,*hdi,*hdj,*hoi,*hoj;
245   PetscInt               *ii,*jj,*iptr,*jptr;
246   PetscInt               cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N;
247   HYPRE_Int              type;
248   PetscErrorCode         ierr;
249 
250   PetscFunctionBegin;
251   comm = PetscObjectComm((PetscObject)A);
252   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type));
253   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
254   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA));
255   M     = hypre_ParCSRMatrixGlobalNumRows(hA);
256   N     = hypre_ParCSRMatrixGlobalNumCols(hA);
257   str   = hypre_ParCSRMatrixFirstRowIndex(hA);
258   stc   = hypre_ParCSRMatrixFirstColDiag(hA);
259   hdiag = hypre_ParCSRMatrixDiag(hA);
260   hoffd = hypre_ParCSRMatrixOffd(hA);
261   dr    = hypre_CSRMatrixNumRows(hdiag);
262   dc    = hypre_CSRMatrixNumCols(hdiag);
263   nnz   = hypre_CSRMatrixNumNonzeros(hdiag);
264   hdi   = hypre_CSRMatrixI(hdiag);
265   hdj   = hypre_CSRMatrixJ(hdiag);
266   hdd   = hypre_CSRMatrixData(hdiag);
267   oc    = hypre_CSRMatrixNumCols(hoffd);
268   nnz  += hypre_CSRMatrixNumNonzeros(hoffd);
269   hoi   = hypre_CSRMatrixI(hoffd);
270   hoj   = hypre_CSRMatrixJ(hoffd);
271   hod   = hypre_CSRMatrixData(hoffd);
272   if (reuse != MAT_REUSE_MATRIX) {
273     PetscInt *aux;
274 
275     /* generate l2g maps for rows and cols */
276     ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr);
277     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
278     ierr = ISDestroy(&is);CHKERRQ(ierr);
279     col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
280     ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr);
281     for (i=0; i<dc; i++) aux[i] = i+stc;
282     for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i];
283     ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
284     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
285     ierr = ISDestroy(&is);CHKERRQ(ierr);
286     /* create MATIS object */
287     ierr = MatCreate(comm,B);CHKERRQ(ierr);
288     ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr);
289     ierr = MatSetType(*B,MATIS);CHKERRQ(ierr);
290     ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr);
291     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
292     ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
293 
294     /* allocate CSR for local matrix */
295     ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr);
296     ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr);
297     ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr);
298   } else {
299     PetscInt  nr;
300     PetscBool done;
301     ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr);
302     ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr);
303     if (nr != dr) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of rows in local mat! %D != %D",nr,dr);
304     if (iptr[nr] < nnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of nonzeros in local mat! reuse %D requested %D",iptr[nr],nnz);
305     ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr);
306   }
307   /* merge local matrices */
308   ii   = iptr;
309   jj   = jptr;
310   aa   = data;
311   *ii  = *(hdi++) + *(hoi++);
312   for (jd=0,jo=0,cum=0; *ii<nnz; cum++) {
313     PetscScalar *aold = aa;
314     PetscInt    *jold = jj,nc = jd+jo;
315     for (; jd<*hdi; jd++) { *jj++ = *hdj++;      *aa++ = *hdd++; }
316     for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
317     *(++ii) = *(hdi++) + *(hoi++);
318     ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr);
319   }
320   for (; cum<dr; cum++) *(++ii) = nnz;
321   if (reuse != MAT_REUSE_MATRIX) {
322     Mat_SeqAIJ* a;
323 
324     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr);
325     ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr);
326     /* hack SeqAIJ */
327     a          = (Mat_SeqAIJ*)(lA->data);
328     a->free_a  = PETSC_TRUE;
329     a->free_ij = PETSC_TRUE;
330     ierr = MatDestroy(&lA);CHKERRQ(ierr);
331   }
332   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
333   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
334   if (reuse == MAT_INPLACE_MATRIX) {
335     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
336   }
337   PetscFunctionReturn(0);
338 }
339 
340 #undef __FUNCT__
341 #define __FUNCT__ "MatConvert_AIJ_HYPRE"
342 PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
343 {
344   Mat_HYPRE      *hB;
345   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
346   PetscErrorCode ierr;
347 
348   PetscFunctionBegin;
349   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
350   if (reuse == MAT_REUSE_MATRIX) {
351     /* always destroy the old matrix and create a new memory;
352        hope this does not churn the memory too much. The problem
353        is I do not know if it is possible to put the matrix back to
354        its initial state so that we can directly copy the values
355        the second time through. */
356     hB = (Mat_HYPRE*)((*B)->data);
357     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
358   } else {
359     ierr = MatCreate(comm,B);CHKERRQ(ierr);
360     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
361     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
362     hB   = (Mat_HYPRE*)((*B)->data);
363   }
364   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
365   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
366   (*B)->preallocated = PETSC_TRUE;
367   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
368   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
369   PetscFunctionReturn(0);
370 }
371 
372 #undef __FUNCT__
373 #define __FUNCT__ "MatConvert_HYPRE_AIJ"
374 static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
375 {
376   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
377   hypre_ParCSRMatrix *parcsr;
378   hypre_CSRMatrix    *hdiag,*hoffd;
379   MPI_Comm           comm;
380   PetscScalar        *da,*oa,*aptr;
381   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
382   PetscInt           i,dnnz,onnz,m,n;
383   HYPRE_Int          type;
384   PetscMPIInt        size;
385   PetscErrorCode     ierr;
386 
387   PetscFunctionBegin;
388   comm = PetscObjectComm((PetscObject)A);
389   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
390   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
391   if (reuse == MAT_REUSE_MATRIX) {
392     PetscBool ismpiaij,isseqaij;
393     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
394     ierr = PetscObjectTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
395     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
396   }
397   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
398 
399   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
400   hdiag = hypre_ParCSRMatrixDiag(parcsr);
401   hoffd = hypre_ParCSRMatrixOffd(parcsr);
402   m     = hypre_CSRMatrixNumRows(hdiag);
403   n     = hypre_CSRMatrixNumCols(hdiag);
404   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
405   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
406   if (reuse == MAT_INITIAL_MATRIX) {
407     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
408     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
409     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
410   } else if (reuse == MAT_REUSE_MATRIX) {
411     PetscInt  nr;
412     PetscBool done;
413     if (size > 1) {
414       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
415 
416       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
417       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);
418       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);
419       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
420     } else {
421       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
422       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
423       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);
424       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
425     }
426   } else { /* MAT_INPLACE_MATRIX */
427     dii = (PetscInt*)hypre_CSRMatrixI(hdiag);
428     djj = (PetscInt*)hypre_CSRMatrixJ(hdiag);
429     da  = hypre_CSRMatrixData(hdiag);
430   }
431   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
432   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
433   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
434   iptr = djj;
435   aptr = da;
436   for (i=0; i<m; i++) {
437     PetscInt nc = dii[i+1]-dii[i];
438     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
439     iptr += nc;
440     aptr += nc;
441   }
442   if (size > 1) {
443     HYPRE_Int *offdj,*coffd;
444 
445     if (reuse == MAT_INITIAL_MATRIX) {
446       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
447       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
448       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
449     } else if (reuse == MAT_REUSE_MATRIX) {
450       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
451       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
452       PetscBool  done;
453 
454       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
455       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);
456       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);
457       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
458     } else { /* MAT_INPLACE_MATRIX */
459       oii = (PetscInt*)hypre_CSRMatrixI(hoffd);
460       ojj = (PetscInt*)hypre_CSRMatrixJ(hoffd);
461       oa  = hypre_CSRMatrixData(hoffd);
462     }
463     ierr  = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
464     offdj = hypre_CSRMatrixJ(hoffd);
465     coffd = hypre_ParCSRMatrixColMapOffd(parcsr);
466     for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]];
467     ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr);
468     iptr = ojj;
469     aptr = oa;
470     for (i=0; i<m; i++) {
471        PetscInt nc = oii[i+1]-oii[i];
472        ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
473        iptr += nc;
474        aptr += nc;
475     }
476     if (reuse == MAT_INITIAL_MATRIX) {
477       Mat_MPIAIJ *b;
478       Mat_SeqAIJ *d,*o;
479 
480       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,B);CHKERRQ(ierr);
481       /* hack MPIAIJ */
482       b          = (Mat_MPIAIJ*)((*B)->data);
483       d          = (Mat_SeqAIJ*)b->A->data;
484       o          = (Mat_SeqAIJ*)b->B->data;
485       d->free_a  = PETSC_TRUE;
486       d->free_ij = PETSC_TRUE;
487       o->free_a  = PETSC_TRUE;
488       o->free_ij = PETSC_TRUE;
489     } else if (reuse == MAT_INPLACE_MATRIX) {
490       Mat T;
491       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,&T);CHKERRQ(ierr);
492       hypre_CSRMatrixI(hdiag)    = NULL;
493       hypre_CSRMatrixJ(hdiag)    = NULL;
494       hypre_CSRMatrixData(hdiag) = NULL;
495       hypre_CSRMatrixI(hoffd)    = NULL;
496       hypre_CSRMatrixJ(hoffd)    = NULL;
497       hypre_CSRMatrixData(hoffd) = NULL;
498       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
499     }
500   } else {
501     oii  = NULL;
502     ojj  = NULL;
503     oa   = NULL;
504     if (reuse == MAT_INITIAL_MATRIX) {
505       Mat_SeqAIJ* b;
506       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
507       /* hack SeqAIJ */
508       b          = (Mat_SeqAIJ*)((*B)->data);
509       b->free_a  = PETSC_TRUE;
510       b->free_ij = PETSC_TRUE;
511     } else if (reuse == MAT_INPLACE_MATRIX) {
512       Mat T;
513       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,&T);CHKERRQ(ierr);
514       hypre_CSRMatrixI(hdiag)    = NULL;
515       hypre_CSRMatrixJ(hdiag)    = NULL;
516       hypre_CSRMatrixData(hdiag) = NULL;
517       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
518     }
519   }
520 
521   /* we have to use hypre_Tfree to free the arrays */
522   if (reuse == MAT_INPLACE_MATRIX) {
523     void *ptrs[6] = {dii,djj,da,oii,ojj,oa};
524     const char *names[6] = {"_hypre_csr_dii",
525                             "_hypre_csr_djj",
526                             "_hypre_csr_da",
527                             "_hypre_csr_oii",
528                             "_hypre_csr_ojj",
529                             "_hypre_csr_oa"};
530     for (i=0; i<6; i++) {
531       PetscContainer c;
532 
533       ierr = PetscContainerCreate(comm,&c);CHKERRQ(ierr);
534       ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr);
535       ierr = PetscContainerSetUserDestroy(c,hypre_array_destroy);CHKERRQ(ierr);
536       ierr = PetscObjectCompose((PetscObject)(*B),names[i],(PetscObject)c);CHKERRQ(ierr);
537       ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
538     }
539   }
540   PetscFunctionReturn(0);
541 }
542 
543 #undef __FUNCT__
544 #define __FUNCT__ "MatAIJGetParCSR_Private"
545 static PetscErrorCode MatAIJGetParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
546 {
547   hypre_ParCSRMatrix *tA;
548   hypre_CSRMatrix    *hdiag,*hoffd;
549   Mat_SeqAIJ         *diag,*offd;
550   PetscInt           *garray,noffd,dnnz,onnz,*row_starts,*col_starts;
551   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
552   PetscBool          ismpiaij,isseqaij;
553   PetscErrorCode     ierr;
554 
555   PetscFunctionBegin;
556   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
557   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
558   if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type);
559   if (ismpiaij) {
560     Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data);
561 
562     diag   = (Mat_SeqAIJ*)a->A->data;
563     offd   = (Mat_SeqAIJ*)a->B->data;
564     garray = a->garray;
565     noffd  = a->B->cmap->N;
566     dnnz   = diag->nz;
567     onnz   = offd->nz;
568   } else {
569     diag   = (Mat_SeqAIJ*)A->data;
570     offd   = NULL;
571     garray = NULL;
572     noffd  = 0;
573     dnnz   = diag->nz;
574     onnz   = 0;
575   }
576 
577   /* create a temporary ParCSR */
578   if (HYPRE_AssumedPartitionCheck()) {
579     PetscMPIInt myid;
580 
581     ierr       = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
582     row_starts = A->rmap->range + myid;
583     col_starts = A->cmap->range + myid;
584   } else {
585     row_starts = A->rmap->range;
586     col_starts = A->cmap->range;
587   }
588   tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,(HYPRE_Int*)row_starts,(HYPRE_Int*)col_starts,noffd,dnnz,onnz);
589   hypre_ParCSRMatrixSetRowStartsOwner(tA,0);
590   hypre_ParCSRMatrixSetColStartsOwner(tA,0);
591 
592   /* set diagonal part */
593   hdiag = hypre_ParCSRMatrixDiag(tA);
594   hypre_CSRMatrixI(hdiag)           = (HYPRE_Int*)diag->i;
595   hypre_CSRMatrixJ(hdiag)           = (HYPRE_Int*)diag->j;
596   hypre_CSRMatrixData(hdiag)        = diag->a;
597   hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz;
598   hypre_CSRMatrixSetRownnz(hdiag);
599   hypre_CSRMatrixSetDataOwner(hdiag,0);
600 
601   /* set offdiagonal part */
602   hoffd = hypre_ParCSRMatrixOffd(tA);
603   if (offd) {
604     hypre_CSRMatrixI(hoffd)           = (HYPRE_Int*)offd->i;
605     hypre_CSRMatrixJ(hoffd)           = (HYPRE_Int*)offd->j;
606     hypre_CSRMatrixData(hoffd)        = offd->a;
607     hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz;
608     hypre_CSRMatrixSetRownnz(hoffd);
609     hypre_CSRMatrixSetDataOwner(hoffd,0);
610     hypre_ParCSRMatrixSetNumNonzeros(tA);
611     hypre_ParCSRMatrixColMapOffd(tA) = (HYPRE_Int*)garray;
612   }
613   *hA = tA;
614   PetscFunctionReturn(0);
615 }
616 
617 #undef __FUNCT__
618 #define __FUNCT__ "MatAIJRestoreParCSR_Private"
619 static PetscErrorCode MatAIJRestoreParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
620 {
621   hypre_CSRMatrix    *hdiag,*hoffd;
622 
623   PetscFunctionBegin;
624   hdiag = hypre_ParCSRMatrixDiag(*hA);
625   hoffd = hypre_ParCSRMatrixOffd(*hA);
626   /* set pointers to NULL before destroying tA */
627   hypre_CSRMatrixI(hdiag)           = NULL;
628   hypre_CSRMatrixJ(hdiag)           = NULL;
629   hypre_CSRMatrixData(hdiag)        = NULL;
630   hypre_CSRMatrixI(hoffd)           = NULL;
631   hypre_CSRMatrixJ(hoffd)           = NULL;
632   hypre_CSRMatrixData(hoffd)        = NULL;
633   hypre_ParCSRMatrixColMapOffd(*hA) = NULL;
634   hypre_ParCSRMatrixDestroy(*hA);
635   *hA = NULL;
636   PetscFunctionReturn(0);
637 }
638 
639 /* calls RAP from BoomerAMG:
640    the resulting ParCSR will not own the column and row starts */
641 #undef __FUNCT__
642 #define __FUNCT__ "MatHYPRE_ParCSR_RAP"
643 static PetscErrorCode MatHYPRE_ParCSR_RAP(hypre_ParCSRMatrix *hR, hypre_ParCSRMatrix *hA,
644                                           hypre_ParCSRMatrix *hP, hypre_ParCSRMatrix **hRAP)
645 {
646   PetscErrorCode ierr;
647   HYPRE_Int      P_owns_col_starts,R_owns_row_starts;
648 
649   PetscFunctionBegin;
650   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(hP);
651   R_owns_row_starts = hypre_ParCSRMatrixOwnsRowStarts(hR);
652   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(hR,hA,hP,hRAP));
653   PetscStackCallStandard(hypre_ParCSRMatrixSetNumNonzeros,(*hRAP));
654   /* hypre_BoomerAMGBuildCoarseOperator steals the col_starts from P and the row_starts from R */
655   hypre_ParCSRMatrixSetRowStartsOwner(*hRAP,0);
656   hypre_ParCSRMatrixSetColStartsOwner(*hRAP,0);
657   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(hP,1);
658   if (R_owns_row_starts) hypre_ParCSRMatrixSetRowStartsOwner(hR,1);
659   PetscFunctionReturn(0);
660 }
661 
662 #undef __FUNCT__
663 #define __FUNCT__ "MatPtAP_AIJ_AIJ_wHYPRE"
664 PETSC_INTERN PetscErrorCode MatPtAP_AIJ_AIJ_wHYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
665 {
666   hypre_ParCSRMatrix *hA,*hP,*hPtAP;
667   PetscErrorCode ierr;
668 
669   PetscFunctionBegin;
670   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
671   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
672   ierr = MatAIJGetParCSR_Private(P,&hP);CHKERRQ(ierr);
673   ierr = MatHYPRE_ParCSR_RAP(hP,hA,hP,&hPtAP);CHKERRQ(ierr);
674   ierr = MatCreateFromParCSR(hPtAP,MATAIJ,PETSC_OWN_POINTER,C);CHKERRQ(ierr);
675   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
676   ierr = MatAIJRestoreParCSR_Private(P,&hP);CHKERRQ(ierr);
677   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
678   PetscFunctionReturn(0);
679 }
680 
681 #undef __FUNCT__
682 #define __FUNCT__ "MatPtAP_AIJ_HYPRE"
683 static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
684 {
685   Mat_HYPRE          *hP = (Mat_HYPRE*)P->data;
686   hypre_ParCSRMatrix *hA = NULL,*Pparcsr,*ptapparcsr;
687   HYPRE_Int          type;
688   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
689   char               mtype[256];
690   PetscErrorCode     ierr;
691 
692   PetscFunctionBegin;
693   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
694   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
695   if (scall == MAT_REUSE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_REUSE_MATRIX");
696   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
697   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
698 
699   /* It looks like we don't need to have the diagonal entries
700      ordered first in the rows of the diagonal part
701      for boomerAMGBuildCoarseOperator to work */
702   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
703   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,hA,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
704   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
705 
706   /* create C depending on mtype */
707   sprintf(mtype,MATAIJ);
708   ierr = PetscOptionsGetString(((PetscObject)A)->options,((PetscObject)A)->prefix,"-matptap_hypre_outtype",mtype,256,NULL);CHKERRQ(ierr);
709   ierr = MatCreateFromParCSR(ptapparcsr,mtype,PETSC_OWN_POINTER,C);CHKERRQ(ierr);
710   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
711   PetscFunctionReturn(0);
712 }
713 
714 #undef __FUNCT__
715 #define __FUNCT__ "MatPtAP_HYPRE_HYPRE"
716 static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
717 {
718   hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr;
719   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data, *hP = (Mat_HYPRE*)P->data;
720   HYPRE_Int          type;
721   PetscErrorCode     ierr;
722 
723   PetscFunctionBegin;
724   if (scall == MAT_REUSE_MATRIX) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Unsupported MAT_REUSE_MATRIX");
725   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
726   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
727   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
728   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
729   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
730   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
731   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
732   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
733   ierr = MatCreateFromParCSR(ptapparcsr,MATHYPRE,PETSC_OWN_POINTER,C);CHKERRQ(ierr);
734   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
735   PetscFunctionReturn(0);
736 }
737 
738 #undef __FUNCT__
739 #define __FUNCT__ "MatMultTranspose_HYPRE"
740 static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
741 {
742   PetscErrorCode ierr;
743 
744   PetscFunctionBegin;
745   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
746   PetscFunctionReturn(0);
747 }
748 
749 #undef __FUNCT__
750 #define __FUNCT__ "MatMult_HYPRE"
751 static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
752 {
753   PetscErrorCode ierr;
754 
755   PetscFunctionBegin;
756   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
757   PetscFunctionReturn(0);
758 }
759 
760 #undef __FUNCT__
761 #define __FUNCT__ "MatHYPRE_MultKernel_Private"
762 static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
763 {
764   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
765   hypre_ParCSRMatrix *parcsr;
766   hypre_ParVector    *hx,*hy;
767   PetscScalar        *ax,*ay,*sax,*say;
768   PetscErrorCode     ierr;
769 
770   PetscFunctionBegin;
771   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
772   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
773   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
774   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
775   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
776   if (trans) {
777     VecHYPRE_ParVectorReplacePointer(hA->x,ay,say);
778     VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax);
779     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
780     VecHYPRE_ParVectorReplacePointer(hA->x,say,ay);
781     VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax);
782   } else {
783     VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax);
784     VecHYPRE_ParVectorReplacePointer(hA->b,ay,say);
785     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
786     VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax);
787     VecHYPRE_ParVectorReplacePointer(hA->b,say,ay);
788   }
789   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
790   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
791   PetscFunctionReturn(0);
792 }
793 
794 #undef __FUNCT__
795 #define __FUNCT__ "MatDestroy_HYPRE"
796 static PetscErrorCode MatDestroy_HYPRE(Mat A)
797 {
798   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
799   PetscErrorCode ierr;
800 
801   PetscFunctionBegin;
802   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
803   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
804   if (hA->ij) {
805     if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL;
806     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
807   }
808   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); }
809   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
810   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr);
811   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr);
812   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr);
813   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPRESetPreallocation_C",NULL);CHKERRQ(ierr);
814   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPREGetParCSR_C",NULL);CHKERRQ(ierr);
815   ierr = PetscFree(A->data);CHKERRQ(ierr);
816   PetscFunctionReturn(0);
817 }
818 
819 #undef __FUNCT__
820 #define __FUNCT__ "MatSetUp_HYPRE"
821 static PetscErrorCode MatSetUp_HYPRE(Mat A)
822 {
823   PetscErrorCode ierr;
824 
825   PetscFunctionBegin;
826   ierr = MatHYPRESetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr);
827   PetscFunctionReturn(0);
828 }
829 
830 #undef __FUNCT__
831 #define __FUNCT__ "MatAssemblyEnd_HYPRE"
832 static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
833 {
834   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
835   Vec                x,b;
836   PetscErrorCode     ierr;
837 
838   PetscFunctionBegin;
839   if (mode == MAT_FLUSH_ASSEMBLY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_FLUSH_ASSEMBLY currently not supported with MATHYPRE");
840   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
841   if (hA->x) PetscFunctionReturn(0);
842   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
843   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
844   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
845   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
846   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
847   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
848   ierr = VecDestroy(&x);CHKERRQ(ierr);
849   ierr = VecDestroy(&b);CHKERRQ(ierr);
850   PetscFunctionReturn(0);
851 }
852 
853 #define MATHYPRE_SCRATCH 2048
854 
855 #undef __FUNCT__
856 #define __FUNCT__ "MatSetValues_HYPRE"
857 PetscErrorCode MatSetValues_HYPRE(Mat A, PetscInt nr, const PetscInt rows[], PetscInt nc, const PetscInt cols[], const PetscScalar v[], InsertMode ins)
858 {
859   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
860   PetscScalar        *vals = (PetscScalar *)v;
861   PetscScalar        sscr[MATHYPRE_SCRATCH];
862   HYPRE_Int          cscr[2][MATHYPRE_SCRATCH];
863   HYPRE_Int          i,nzc;
864   PetscErrorCode     ierr;
865 
866   PetscFunctionBegin;
867   for (i=0,nzc=0;i<nc;i++) {
868     if (cols[i] >= 0) {
869       cscr[0][nzc  ] = cols[i];
870       cscr[1][nzc++] = i;
871     }
872   }
873   if (!nzc) PetscFunctionReturn(0);
874 
875   if (ins == ADD_VALUES) {
876     for (i=0;i<nr;i++) {
877       if (rows[i] >= 0) {
878         PetscInt j;
879         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
880         PetscStackCallStandard(HYPRE_IJMatrixAddToValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
881       }
882       vals += nc;
883     }
884   } else { /* INSERT_VALUES */
885 #if defined(PETSC_USE_DEBUG)
886     /* Insert values cannot be used to insert offproc entries */
887     PetscInt rst,ren;
888     ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
889     for (i=0;i<nr;i++)
890       if (rows[i] >= 0 && (rows[i] < rst || rows[i] >= ren)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot use INSERT_VALUES for off-proc entries with MatHYPRE. Use ADD_VALUES instead");
891 #endif
892     for (i=0;i<nr;i++) {
893       if (rows[i] >= 0) {
894         PetscInt j;
895         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
896         PetscStackCallStandard(HYPRE_IJMatrixSetValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
897       }
898       vals += nc;
899     }
900   }
901   PetscFunctionReturn(0);
902 }
903 
904 #undef __FUNCT__
905 #define __FUNCT__ "MatHYPRESetPreallocation_HYPRE"
906 static PetscErrorCode MatHYPRESetPreallocation_HYPRE(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
907 {
908   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
909   HYPRE_Int          *hdnnz,*honnz;
910   PetscInt           i,rs,re,cs,ce,bs;
911   PetscMPIInt        size;
912   PetscErrorCode     ierr;
913 
914   PetscFunctionBegin;
915   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
916   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
917   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
918   rs   = A->rmap->rstart;
919   re   = A->rmap->rend;
920   cs   = A->cmap->rstart;
921   ce   = A->cmap->rend;
922   if (!hA->ij) {
923     PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rs,re-1,cs,ce-1,&hA->ij));
924     PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
925   } else {
926     HYPRE_Int hrs,hre,hcs,hce;
927     PetscStackCallStandard(HYPRE_IJMatrixGetLocalRange,(hA->ij,&hrs,&hre,&hcs,&hce));
928     if (hre-hrs+1 != re -rs) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent local rows: IJMatrix [%D,%D), PETSc [%D,%d)",hrs,hre+1,rs,re);
929     if (hce-hcs+1 != ce -cs) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inconsistent local cols: IJMatrix [%D,%D), PETSc [%D,%d)",hcs,hce+1,cs,ce);
930   }
931   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
932 
933   if (dnz == PETSC_DEFAULT || dnz == PETSC_DECIDE) dnz = 10*bs;
934   if (onz == PETSC_DEFAULT || onz == PETSC_DECIDE) onz = 10*bs;
935 
936   if (!dnnz) {
937     ierr = PetscMalloc1(A->rmap->n,&hdnnz);CHKERRQ(ierr);
938     for (i=0;i<A->rmap->n;i++) hdnnz[i] = dnz;
939   } else {
940     hdnnz = (HYPRE_Int*)dnnz;
941   }
942   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
943   if (size > 1) {
944     if (!onnz) {
945       ierr = PetscMalloc1(A->rmap->n,&honnz);CHKERRQ(ierr);
946       for (i=0;i<A->rmap->n;i++) honnz[i] = onz;
947     } else {
948       honnz = (HYPRE_Int*)onnz;
949     }
950     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(hA->ij,hdnnz,honnz));
951   } else {
952     honnz = NULL;
953     PetscStackCallStandard(HYPRE_IJMatrixSetRowSizes,(hA->ij,hdnnz));
954   }
955   if (!dnnz) {
956     ierr = PetscFree(hdnnz);CHKERRQ(ierr);
957   }
958   if (!onnz && honnz) {
959     ierr = PetscFree(honnz);CHKERRQ(ierr);
960   }
961   A->preallocated = PETSC_TRUE;
962 
963   /* SetDiagOffdSizes sets hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0 */
964   {
965     hypre_AuxParCSRMatrix *aux_matrix;
966     aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij);
967     hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1;
968   }
969   PetscFunctionReturn(0);
970 }
971 
972 /*@C
973    MatHYPRESetPreallocation - Preallocates memory for a sparse parallel matrix in HYPRE IJ format
974 
975    Collective on Mat
976 
977    Input Parameters:
978 +  A - the matrix
979 .  dnz  - number of nonzeros per row in DIAGONAL portion of local submatrix
980           (same value is used for all local rows)
981 .  dnnz - array containing the number of nonzeros in the various rows of the
982           DIAGONAL portion of the local submatrix (possibly different for each row)
983           or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure.
984           The size of this array is equal to the number of local rows, i.e 'm'.
985           For matrices that will be factored, you must leave room for (and set)
986           the diagonal entry even if it is zero.
987 .  onz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
988           submatrix (same value is used for all local rows).
989 -  onnz - array containing the number of nonzeros in the various rows of the
990           OFF-DIAGONAL portion of the local submatrix (possibly different for
991           each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero
992           structure. The size of this array is equal to the number
993           of local rows, i.e 'm'.
994 
995    Notes: If the *nnz parameter is given then the *nz parameter is ignored; for sequential matrices, onz and onnz are ignored.
996 
997    Level: intermediate
998 
999 .keywords: matrix, aij, compressed row, sparse, parallel
1000 
1001 .seealso: MatCreate(), MatMPIAIJSetPreallocation, MATHYPRE
1002 @*/
1003 #undef __FUNCT__
1004 #define __FUNCT__ "MatHYPRESetPreallocation"
1005 PetscErrorCode MatHYPRESetPreallocation(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1006 {
1007   PetscErrorCode ierr;
1008 
1009   PetscFunctionBegin;
1010   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1011   PetscValidType(A,1);
1012   ierr = PetscTryMethod(A,"MatHYPRESetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(A,dnz,dnnz,onz,onnz));CHKERRQ(ierr);
1013   PetscFunctionReturn(0);
1014 }
1015 
1016 /*
1017    MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix
1018 
1019    Collective
1020 
1021    Input Parameters:
1022 +  vparcsr  - the pointer to the hypre_ParCSRMatrix
1023 .  mtype    - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported.
1024 -  copymode - PETSc copying options
1025 
1026    Output Parameter:
1027 .  A  - the matrix
1028 
1029    Level: intermediate
1030 
1031 .seealso: MatHYPRE, PetscCopyMode
1032 */
1033 #undef __FUNCT__
1034 #define __FUNCT__ "MatCreateFromParCSR"
1035 PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(hypre_ParCSRMatrix *vparcsr, MatType mtype, PetscCopyMode copymode, Mat* A)
1036 {
1037   Mat                   T;
1038   Mat_HYPRE             *hA;
1039   hypre_ParCSRMatrix    *parcsr;
1040   MPI_Comm              comm;
1041   PetscInt              rstart,rend,cstart,cend,M,N;
1042   PetscBool             isseqaij,ismpiaij,isaij,ishyp,isis;
1043   PetscErrorCode        ierr;
1044 
1045   PetscFunctionBegin;
1046   parcsr = (hypre_ParCSRMatrix *)vparcsr;
1047   comm   = hypre_ParCSRMatrixComm(parcsr);
1048   ierr   = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
1049   ierr   = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
1050   ierr   = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr);
1051   ierr   = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr);
1052   ierr   = PetscStrcmp(mtype,MATIS,&isis);CHKERRQ(ierr);
1053   isaij  = (PetscBool)(isseqaij || ismpiaij || isaij);
1054   if (!isaij && !ishyp && !isis) SETERRQ6(comm,PETSC_ERR_SUP,"Unsupported MatType %s! Supported types are %s, %s, %s, %s, and %s",mtype,MATAIJ,MATSEQAIJ,MATMPIAIJ,MATIS,MATHYPRE);
1055   if (ishyp && copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES");
1056 
1057   /* access ParCSRMatrix */
1058   rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr);
1059   rend   = hypre_ParCSRMatrixLastRowIndex(parcsr);
1060   cstart = hypre_ParCSRMatrixFirstColDiag(parcsr);
1061   cend   = hypre_ParCSRMatrixLastColDiag(parcsr);
1062   M      = hypre_ParCSRMatrixGlobalNumRows(parcsr);
1063   N      = hypre_ParCSRMatrixGlobalNumCols(parcsr);
1064 
1065   /* fix for empty local rows/columns */
1066   if (rend < rstart) rend = rstart;
1067   if (cend < cstart) cend = cstart;
1068 
1069   /* create PETSc matrix with MatHYPRE */
1070   ierr = MatCreate(comm,&T);CHKERRQ(ierr);
1071   ierr = MatSetSizes(T,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr);
1072   ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr);
1073   hA   = (Mat_HYPRE*)(T->data);
1074 
1075   /* create HYPRE_IJMatrix */
1076   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
1077 
1078   /* set ParCSR object */
1079   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1080   hypre_IJMatrixObject(hA->ij) = parcsr;
1081   T->preallocated = PETSC_TRUE;
1082 
1083   /* set assembled flag */
1084   hypre_IJMatrixAssembleFlag(hA->ij) = 1;
1085   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1086   if (ishyp) {
1087     PetscMPIInt myid = 0;
1088 
1089     /* make sure we always have row_starts and col_starts available */
1090     if (HYPRE_AssumedPartitionCheck()) {
1091       ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
1092     }
1093     if (!hypre_ParCSRMatrixOwnsColStarts(parcsr)) {
1094       PetscLayout map;
1095 
1096       ierr = MatGetLayouts(T,NULL,&map);CHKERRQ(ierr);
1097       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1098       hypre_ParCSRMatrixColStarts(parcsr) = map->range + myid;
1099     }
1100     if (!hypre_ParCSRMatrixOwnsRowStarts(parcsr)) {
1101       PetscLayout map;
1102 
1103       ierr = MatGetLayouts(T,&map,NULL);CHKERRQ(ierr);
1104       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1105       hypre_ParCSRMatrixRowStarts(parcsr) = map->range + myid;
1106     }
1107     /* prevent from freeing the pointer */
1108     if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE;
1109     *A = T;
1110     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1111     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1112   } else if (isaij) {
1113     if (copymode != PETSC_OWN_POINTER) {
1114       /* prevent from freeing the pointer */
1115       hA->inner_free = PETSC_FALSE;
1116       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
1117       ierr = MatDestroy(&T);CHKERRQ(ierr);
1118     } else { /* AIJ return type with PETSC_OWN_POINTER */
1119       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1120       *A   = T;
1121     }
1122   } else if (isis) {
1123     ierr = MatConvert_HYPRE_IS(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
1124     if (copymode != PETSC_OWN_POINTER) hA->inner_free = PETSC_FALSE;
1125     ierr = MatDestroy(&T);CHKERRQ(ierr);
1126   }
1127   PetscFunctionReturn(0);
1128 }
1129 
1130 #undef __FUNCT__
1131 #define __FUNCT__ "MatHYPREGetParCSR_HYPRE"
1132 PetscErrorCode MatHYPREGetParCSR_HYPRE(Mat A, hypre_ParCSRMatrix **parcsr)
1133 {
1134   Mat_HYPRE*            hA = (Mat_HYPRE*)A->data;
1135   HYPRE_Int             type;
1136   PetscErrorCode        ierr;
1137 
1138   PetscFunctionBegin;
1139   if (!hA->ij) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"HYPRE_IJMatrix not present");
1140   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
1141   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"HYPRE_IJMatrix is not of type HYPRE_PARCSR");
1142   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)parcsr));
1143   PetscFunctionReturn(0);
1144 }
1145 
1146 /*
1147    MatHYPREGetParCSR - Gets the pointer to the ParCSR matrix
1148 
1149    Not collective
1150 
1151    Input Parameters:
1152 +  A  - the MATHYPRE object
1153 
1154    Output Parameter:
1155 .  parcsr  - the pointer to the hypre_ParCSRMatrix
1156 
1157    Level: intermediate
1158 
1159 .seealso: MatHYPRE, PetscCopyMode
1160 */
1161 #undef __FUNCT__
1162 #define __FUNCT__ "MatHYPREGetParCSR"
1163 PetscErrorCode MatHYPREGetParCSR(Mat A, hypre_ParCSRMatrix **parcsr)
1164 {
1165   PetscErrorCode ierr;
1166 
1167   PetscFunctionBegin;
1168   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1169   PetscValidType(A,1);
1170   ierr = PetscUseMethod(A,"MatHYPREGetParCSR_C",(Mat,hypre_ParCSRMatrix**),(A,parcsr));CHKERRQ(ierr);
1171   PetscFunctionReturn(0);
1172 }
1173 
1174 #undef __FUNCT__
1175 #define __FUNCT__ "MatCreate_HYPRE"
1176 PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
1177 {
1178   Mat_HYPRE      *hB;
1179   PetscErrorCode ierr;
1180 
1181   PetscFunctionBegin;
1182   ierr           = PetscNewLog(B,&hB);CHKERRQ(ierr);
1183   hB->inner_free = PETSC_TRUE;
1184 
1185   B->data       = (void*)hB;
1186   B->rmap->bs   = 1;
1187   B->assembled  = PETSC_FALSE;
1188 
1189   B->ops->mult          = MatMult_HYPRE;
1190   B->ops->multtranspose = MatMultTranspose_HYPRE;
1191   B->ops->setup         = MatSetUp_HYPRE;
1192   B->ops->destroy       = MatDestroy_HYPRE;
1193   B->ops->assemblyend   = MatAssemblyEnd_HYPRE;
1194   B->ops->ptap          = MatPtAP_HYPRE_HYPRE;
1195   B->ops->setvalues     = MatSetValues_HYPRE;
1196 
1197   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
1198   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
1199   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
1200   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr);
1201   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1202   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1203   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPRESetPreallocation_C",MatHYPRESetPreallocation_HYPRE);CHKERRQ(ierr);
1204   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPREGetParCSR_C",MatHYPREGetParCSR_HYPRE);CHKERRQ(ierr);
1205   PetscFunctionReturn(0);
1206 }
1207 
1208 #undef __FUNCT__
1209 #define __FUNCT__ "hypre_array_destroy"
1210 static PetscErrorCode hypre_array_destroy(void *ptr)
1211 {
1212    PetscFunctionBegin;
1213    hypre_TFree(ptr);
1214    PetscFunctionReturn(0);
1215 }
1216 
1217 #if 0
1218 /*
1219     Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format
1220 
1221     This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first
1222     which will corrupt the PETSc data structure if we did this. Need a work around to this problem.
1223 */
1224 #include <_hypre_IJ_mv.h>
1225 #include <HYPRE_IJ_mv.h>
1226 
1227 #undef __FUNCT__
1228 #define __FUNCT__ "MatHYPRE_IJMatrixLink"
1229 PetscErrorCode MatHYPRE_IJMatrixLink(Mat A, HYPRE_IJMatrix *ij)
1230 {
1231   PetscErrorCode        ierr;
1232   PetscInt              rstart,rend,cstart,cend;
1233   PetscBool             flg;
1234   hypre_AuxParCSRMatrix *aux_matrix;
1235 
1236   PetscFunctionBegin;
1237   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1238   PetscValidType(A,1);
1239   PetscValidPointer(ij,2);
1240   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
1241   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices");
1242   ierr = MatSetUp(A);CHKERRQ(ierr);
1243 
1244   ierr   = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
1245   rstart = A->rmap->rstart;
1246   rend   = A->rmap->rend;
1247   cstart = A->cmap->rstart;
1248   cend   = A->cmap->rend;
1249   PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij));
1250   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
1251 
1252   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij));
1253   PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij));
1254 
1255   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
1256 
1257   /* this is the Hack part where we monkey directly with the hypre datastructures */
1258 
1259   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij));
1260   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
1261   PetscFunctionReturn(0);
1262 }
1263 #endif
1264