xref: /petsc/src/mat/impls/hypre/mhypre.c (revision c69f721f3e436f4ba820575611f0ca1fcc5fe181)
1 
2 /*
3     Creates hypre ijmatrix from PETSc matrix
4 */
5 
6 #include <petscmathypre.h>
7 #include <petsc/private/matimpl.h>
8 #include <../src/mat/impls/hypre/mhypre.h>
9 #include <../src/mat/impls/aij/mpi/mpiaij.h>
10 #include <../src/vec/vec/impls/hypre/vhyp.h>
11 #include <HYPRE.h>
12 #include <HYPRE_utilities.h>
13 #include <_hypre_parcsr_ls.h>
14 #include <_hypre_sstruct_ls.h>
15 
16 static PetscErrorCode MatHYPRE_CreateFromMat(Mat,Mat_HYPRE*);
17 static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat,Mat,HYPRE_IJMatrix);
18 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
19 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
20 static PetscErrorCode MatHYPRE_MultKernel_Private(Mat,Vec,Vec,PetscBool);
21 static PetscErrorCode hypre_array_destroy(void*);
22 PetscErrorCode MatSetValues_HYPRE(Mat, PetscInt,const PetscInt[],PetscInt,const PetscInt[],const PetscScalar[],InsertMode ins);
23 
24 static PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o, HYPRE_IJMatrix ij)
25 {
26   PetscErrorCode ierr;
27   PetscInt       i,n_d,n_o;
28   const PetscInt *ia_d,*ia_o;
29   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
30   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
31 
32   PetscFunctionBegin;
33   if (A_d) { /* determine number of nonzero entries in local diagonal part */
34     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
35     if (done_d) {
36       ierr = PetscMalloc1(n_d,&nnz_d);CHKERRQ(ierr);
37       for (i=0; i<n_d; i++) {
38         nnz_d[i] = ia_d[i+1] - ia_d[i];
39       }
40     }
41     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
42   }
43   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
44     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
45     if (done_o) {
46       ierr = PetscMalloc1(n_o,&nnz_o);CHKERRQ(ierr);
47       for (i=0; i<n_o; i++) {
48         nnz_o[i] = ia_o[i+1] - ia_o[i];
49       }
50     }
51     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
52   }
53   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
54     if (!done_o) { /* only diagonal part */
55       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
56       for (i=0; i<n_d; i++) {
57         nnz_o[i] = 0;
58       }
59     }
60     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)nnz_o));
61     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
62     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
63   }
64   PetscFunctionReturn(0);
65 }
66 
67 static PetscErrorCode MatHYPRE_CreateFromMat(Mat A, Mat_HYPRE *hA)
68 {
69   PetscErrorCode ierr;
70   PetscInt       rstart,rend,cstart,cend;
71 
72   PetscFunctionBegin;
73   ierr   = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
74   ierr   = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
75   rstart = A->rmap->rstart;
76   rend   = A->rmap->rend;
77   cstart = A->cmap->rstart;
78   cend   = A->cmap->rend;
79   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
80   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
81   {
82     PetscBool      same;
83     Mat            A_d,A_o;
84     const PetscInt *colmap;
85     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
86     if (same) {
87       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
88       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
89       PetscFunctionReturn(0);
90     }
91     ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
92     if (same) {
93       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
94       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,hA->ij);CHKERRQ(ierr);
95       PetscFunctionReturn(0);
96     }
97     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
98     if (same) {
99       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
100       PetscFunctionReturn(0);
101     }
102     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
103     if (same) {
104       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,hA->ij);CHKERRQ(ierr);
105       PetscFunctionReturn(0);
106     }
107   }
108   PetscFunctionReturn(0);
109 }
110 
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     if (ncols) {
137       PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
138     }
139     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
140   }
141   PetscFunctionReturn(0);
142 }
143 
144 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A, HYPRE_IJMatrix ij)
145 {
146   PetscErrorCode        ierr;
147   Mat_SeqAIJ            *pdiag = (Mat_SeqAIJ*)A->data;
148   HYPRE_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 static PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A, HYPRE_IJMatrix ij)
172 {
173   PetscErrorCode        ierr;
174   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
175   Mat_SeqAIJ            *pdiag,*poffd;
176   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
177   HYPRE_Int             type;
178   hypre_ParCSRMatrix    *par_matrix;
179   hypre_AuxParCSRMatrix *aux_matrix;
180   hypre_CSRMatrix       *hdiag,*hoffd;
181 
182   PetscFunctionBegin;
183   pdiag = (Mat_SeqAIJ*) pA->A->data;
184   poffd = (Mat_SeqAIJ*) pA->B->data;
185   /* cstart is only valid for square MPIAIJ layed out in the usual way */
186   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
187 
188   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
189   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(ij,&type));
190   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
191   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(ij,(void**)&par_matrix));
192   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
193   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
194 
195   /*
196        this is the Hack part where we monkey directly with the hypre datastructures
197   */
198   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
199   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
200   jj  = (PetscInt*)hdiag->j;
201   pjj = (PetscInt*)pdiag->j;
202   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
203   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
204   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
205   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
206      If we hacked a hypre a bit more we might be able to avoid this step */
207   jj  = (PetscInt*) hoffd->j;
208   pjj = (PetscInt*) poffd->j;
209   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
210   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
211 
212   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
213   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
214   PetscFunctionReturn(0);
215 }
216 
217 static PetscErrorCode MatConvert_HYPRE_IS(Mat A, MatType mtype, MatReuse reuse, Mat* B)
218 {
219   Mat_HYPRE*             mhA = (Mat_HYPRE*)(A->data);
220   Mat                    lA;
221   ISLocalToGlobalMapping rl2g,cl2g;
222   IS                     is;
223   hypre_ParCSRMatrix     *hA;
224   hypre_CSRMatrix        *hdiag,*hoffd;
225   MPI_Comm               comm;
226   PetscScalar            *hdd,*hod,*aa,*data;
227   HYPRE_Int              *col_map_offd,*hdi,*hdj,*hoi,*hoj;
228   PetscInt               *ii,*jj,*iptr,*jptr;
229   PetscInt               cum,dr,dc,oc,str,stc,nnz,i,jd,jo,M,N;
230   HYPRE_Int              type;
231   PetscErrorCode         ierr;
232 
233   PetscFunctionBegin;
234   comm = PetscObjectComm((PetscObject)A);
235   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(mhA->ij,&type));
236   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
237   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(mhA->ij,(void**)&hA));
238   M     = hypre_ParCSRMatrixGlobalNumRows(hA);
239   N     = hypre_ParCSRMatrixGlobalNumCols(hA);
240   str   = hypre_ParCSRMatrixFirstRowIndex(hA);
241   stc   = hypre_ParCSRMatrixFirstColDiag(hA);
242   hdiag = hypre_ParCSRMatrixDiag(hA);
243   hoffd = hypre_ParCSRMatrixOffd(hA);
244   dr    = hypre_CSRMatrixNumRows(hdiag);
245   dc    = hypre_CSRMatrixNumCols(hdiag);
246   nnz   = hypre_CSRMatrixNumNonzeros(hdiag);
247   hdi   = hypre_CSRMatrixI(hdiag);
248   hdj   = hypre_CSRMatrixJ(hdiag);
249   hdd   = hypre_CSRMatrixData(hdiag);
250   oc    = hypre_CSRMatrixNumCols(hoffd);
251   nnz  += hypre_CSRMatrixNumNonzeros(hoffd);
252   hoi   = hypre_CSRMatrixI(hoffd);
253   hoj   = hypre_CSRMatrixJ(hoffd);
254   hod   = hypre_CSRMatrixData(hoffd);
255   if (reuse != MAT_REUSE_MATRIX) {
256     PetscInt *aux;
257 
258     /* generate l2g maps for rows and cols */
259     ierr = ISCreateStride(comm,dr,str,1,&is);CHKERRQ(ierr);
260     ierr = ISLocalToGlobalMappingCreateIS(is,&rl2g);CHKERRQ(ierr);
261     ierr = ISDestroy(&is);CHKERRQ(ierr);
262     col_map_offd = hypre_ParCSRMatrixColMapOffd(hA);
263     ierr = PetscMalloc1(dc+oc,&aux);CHKERRQ(ierr);
264     for (i=0; i<dc; i++) aux[i] = i+stc;
265     for (i=0; i<oc; i++) aux[i+dc] = col_map_offd[i];
266     ierr = ISCreateGeneral(comm,dc+oc,aux,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
267     ierr = ISLocalToGlobalMappingCreateIS(is,&cl2g);CHKERRQ(ierr);
268     ierr = ISDestroy(&is);CHKERRQ(ierr);
269     /* create MATIS object */
270     ierr = MatCreate(comm,B);CHKERRQ(ierr);
271     ierr = MatSetSizes(*B,dr,dc,M,N);CHKERRQ(ierr);
272     ierr = MatSetType(*B,MATIS);CHKERRQ(ierr);
273     ierr = MatSetLocalToGlobalMapping(*B,rl2g,cl2g);CHKERRQ(ierr);
274     ierr = ISLocalToGlobalMappingDestroy(&rl2g);CHKERRQ(ierr);
275     ierr = ISLocalToGlobalMappingDestroy(&cl2g);CHKERRQ(ierr);
276 
277     /* allocate CSR for local matrix */
278     ierr = PetscMalloc1(dr+1,&iptr);CHKERRQ(ierr);
279     ierr = PetscMalloc1(nnz,&jptr);CHKERRQ(ierr);
280     ierr = PetscMalloc1(nnz,&data);CHKERRQ(ierr);
281   } else {
282     PetscInt  nr;
283     PetscBool done;
284     ierr = MatISGetLocalMat(*B,&lA);CHKERRQ(ierr);
285     ierr = MatGetRowIJ(lA,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&iptr,(const PetscInt**)&jptr,&done);CHKERRQ(ierr);
286     if (nr != dr) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of rows in local mat! %D != %D",nr,dr);
287     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);
288     ierr = MatSeqAIJGetArray(lA,&data);CHKERRQ(ierr);
289   }
290   /* merge local matrices */
291   ii   = iptr;
292   jj   = jptr;
293   aa   = data;
294   *ii  = *(hdi++) + *(hoi++);
295   for (jd=0,jo=0,cum=0; *ii<nnz; cum++) {
296     PetscScalar *aold = aa;
297     PetscInt    *jold = jj,nc = jd+jo;
298     for (; jd<*hdi; jd++) { *jj++ = *hdj++;      *aa++ = *hdd++; }
299     for (; jo<*hoi; jo++) { *jj++ = *hoj++ + dc; *aa++ = *hod++; }
300     *(++ii) = *(hdi++) + *(hoi++);
301     ierr = PetscSortIntWithScalarArray(jd+jo-nc,jold,aold);CHKERRQ(ierr);
302   }
303   for (; cum<dr; cum++) *(++ii) = nnz;
304   if (reuse != MAT_REUSE_MATRIX) {
305     Mat_SeqAIJ* a;
306 
307     ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,dr,dc+oc,iptr,jptr,data,&lA);CHKERRQ(ierr);
308     ierr = MatISSetLocalMat(*B,lA);CHKERRQ(ierr);
309     /* hack SeqAIJ */
310     a          = (Mat_SeqAIJ*)(lA->data);
311     a->free_a  = PETSC_TRUE;
312     a->free_ij = PETSC_TRUE;
313     ierr = MatDestroy(&lA);CHKERRQ(ierr);
314   }
315   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
316   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
317   if (reuse == MAT_INPLACE_MATRIX) {
318     ierr = MatHeaderReplace(A,B);CHKERRQ(ierr);
319   }
320   PetscFunctionReturn(0);
321 }
322 
323 PETSC_INTERN PetscErrorCode MatConvert_AIJ_HYPRE(Mat A, MatType type, MatReuse reuse, Mat *B)
324 {
325   Mat_HYPRE      *hB;
326   MPI_Comm       comm = PetscObjectComm((PetscObject)A);
327   PetscErrorCode ierr;
328 
329   PetscFunctionBegin;
330   if (reuse == MAT_INPLACE_MATRIX) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported MAT_INPLACE_MATRIX");
331   if (reuse == MAT_REUSE_MATRIX) {
332     /* always destroy the old matrix and create a new memory;
333        hope this does not churn the memory too much. The problem
334        is I do not know if it is possible to put the matrix back to
335        its initial state so that we can directly copy the values
336        the second time through. */
337     hB = (Mat_HYPRE*)((*B)->data);
338     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hB->ij));
339   } else {
340     ierr = MatCreate(comm,B);CHKERRQ(ierr);
341     ierr = MatSetType(*B,MATHYPRE);CHKERRQ(ierr);
342     ierr = MatSetSizes(*B,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
343     hB   = (Mat_HYPRE*)((*B)->data);
344   }
345   ierr = MatHYPRE_CreateFromMat(A,hB);CHKERRQ(ierr);
346   ierr = MatHYPRE_IJMatrixCopy(A,hB->ij);CHKERRQ(ierr);
347   (*B)->preallocated = PETSC_TRUE;
348   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
349   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
350   PetscFunctionReturn(0);
351 }
352 
353 static PetscErrorCode MatConvert_HYPRE_AIJ(Mat A, MatType mtype, MatReuse reuse, Mat *B)
354 {
355   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
356   hypre_ParCSRMatrix *parcsr;
357   hypre_CSRMatrix    *hdiag,*hoffd;
358   MPI_Comm           comm;
359   PetscScalar        *da,*oa,*aptr;
360   PetscInt           *dii,*djj,*oii,*ojj,*iptr;
361   PetscInt           i,dnnz,onnz,m,n;
362   HYPRE_Int          type;
363   PetscMPIInt        size;
364   PetscErrorCode     ierr;
365 
366   PetscFunctionBegin;
367   comm = PetscObjectComm((PetscObject)A);
368   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
369   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
370   if (reuse == MAT_REUSE_MATRIX) {
371     PetscBool ismpiaij,isseqaij;
372     ierr = PetscObjectTypeCompare((PetscObject)*B,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
373     ierr = PetscObjectBaseTypeCompare((PetscObject)*B,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
374     if (!ismpiaij && !isseqaij) SETERRQ(comm,PETSC_ERR_SUP,"Only MATMPIAIJ or MATSEQAIJ are supported");
375   }
376   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
377 
378   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
379   hdiag = hypre_ParCSRMatrixDiag(parcsr);
380   hoffd = hypre_ParCSRMatrixOffd(parcsr);
381   m     = hypre_CSRMatrixNumRows(hdiag);
382   n     = hypre_CSRMatrixNumCols(hdiag);
383   dnnz  = hypre_CSRMatrixNumNonzeros(hdiag);
384   onnz  = hypre_CSRMatrixNumNonzeros(hoffd);
385   if (reuse == MAT_INITIAL_MATRIX) {
386     ierr = PetscMalloc1(m+1,&dii);CHKERRQ(ierr);
387     ierr = PetscMalloc1(dnnz,&djj);CHKERRQ(ierr);
388     ierr = PetscMalloc1(dnnz,&da);CHKERRQ(ierr);
389   } else if (reuse == MAT_REUSE_MATRIX) {
390     PetscInt  nr;
391     PetscBool done;
392     if (size > 1) {
393       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
394 
395       ierr = MatGetRowIJ(b->A,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
396       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);
397       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);
398       ierr = MatSeqAIJGetArray(b->A,&da);CHKERRQ(ierr);
399     } else {
400       ierr = MatGetRowIJ(*B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&dii,(const PetscInt**)&djj,&done);CHKERRQ(ierr);
401       if (nr != m) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_USER,"Cannot reuse mat: invalid number of local rows! %D != %D",nr,m);
402       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);
403       ierr = MatSeqAIJGetArray(*B,&da);CHKERRQ(ierr);
404     }
405   } else { /* MAT_INPLACE_MATRIX */
406     dii = (PetscInt*)hypre_CSRMatrixI(hdiag);
407     djj = (PetscInt*)hypre_CSRMatrixJ(hdiag);
408     da  = hypre_CSRMatrixData(hdiag);
409   }
410   ierr = PetscMemcpy(dii,hypre_CSRMatrixI(hdiag),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
411   ierr = PetscMemcpy(djj,hypre_CSRMatrixJ(hdiag),dnnz*sizeof(PetscInt));CHKERRQ(ierr);
412   ierr = PetscMemcpy(da,hypre_CSRMatrixData(hdiag),dnnz*sizeof(PetscScalar));CHKERRQ(ierr);
413   iptr = djj;
414   aptr = da;
415   for (i=0; i<m; i++) {
416     PetscInt nc = dii[i+1]-dii[i];
417     ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
418     iptr += nc;
419     aptr += nc;
420   }
421   if (size > 1) {
422     HYPRE_Int *offdj,*coffd;
423 
424     if (reuse == MAT_INITIAL_MATRIX) {
425       ierr = PetscMalloc1(m+1,&oii);CHKERRQ(ierr);
426       ierr = PetscMalloc1(onnz,&ojj);CHKERRQ(ierr);
427       ierr = PetscMalloc1(onnz,&oa);CHKERRQ(ierr);
428     } else if (reuse == MAT_REUSE_MATRIX) {
429       Mat_MPIAIJ *b = (Mat_MPIAIJ*)((*B)->data);
430       PetscInt   nr,hr = hypre_CSRMatrixNumRows(hoffd);
431       PetscBool  done;
432 
433       ierr = MatGetRowIJ(b->B,0,PETSC_FALSE,PETSC_FALSE,&nr,(const PetscInt**)&oii,(const PetscInt**)&ojj,&done);CHKERRQ(ierr);
434       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);
435       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);
436       ierr = MatSeqAIJGetArray(b->B,&oa);CHKERRQ(ierr);
437     } else { /* MAT_INPLACE_MATRIX */
438       oii = (PetscInt*)hypre_CSRMatrixI(hoffd);
439       ojj = (PetscInt*)hypre_CSRMatrixJ(hoffd);
440       oa  = hypre_CSRMatrixData(hoffd);
441     }
442     ierr  = PetscMemcpy(oii,hypre_CSRMatrixI(hoffd),(m+1)*sizeof(PetscInt));CHKERRQ(ierr);
443     offdj = hypre_CSRMatrixJ(hoffd);
444     coffd = hypre_ParCSRMatrixColMapOffd(parcsr);
445     for (i=0; i<onnz; i++) ojj[i] = coffd[offdj[i]];
446     ierr = PetscMemcpy(oa,hypre_CSRMatrixData(hoffd),onnz*sizeof(PetscScalar));CHKERRQ(ierr);
447     iptr = ojj;
448     aptr = oa;
449     for (i=0; i<m; i++) {
450        PetscInt nc = oii[i+1]-oii[i];
451        ierr = PetscSortIntWithScalarArray(nc,iptr,aptr);CHKERRQ(ierr);
452        iptr += nc;
453        aptr += nc;
454     }
455     if (reuse == MAT_INITIAL_MATRIX) {
456       Mat_MPIAIJ *b;
457       Mat_SeqAIJ *d,*o;
458 
459       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,B);CHKERRQ(ierr);
460       /* hack MPIAIJ */
461       b          = (Mat_MPIAIJ*)((*B)->data);
462       d          = (Mat_SeqAIJ*)b->A->data;
463       o          = (Mat_SeqAIJ*)b->B->data;
464       d->free_a  = PETSC_TRUE;
465       d->free_ij = PETSC_TRUE;
466       o->free_a  = PETSC_TRUE;
467       o->free_ij = PETSC_TRUE;
468     } else if (reuse == MAT_INPLACE_MATRIX) {
469       Mat T;
470       ierr = MatCreateMPIAIJWithSplitArrays(comm,m,n,PETSC_DECIDE,PETSC_DECIDE,dii,djj,da,oii,ojj,oa,&T);CHKERRQ(ierr);
471       hypre_CSRMatrixI(hdiag)    = NULL;
472       hypre_CSRMatrixJ(hdiag)    = NULL;
473       hypre_CSRMatrixData(hdiag) = NULL;
474       hypre_CSRMatrixI(hoffd)    = NULL;
475       hypre_CSRMatrixJ(hoffd)    = NULL;
476       hypre_CSRMatrixData(hoffd) = NULL;
477       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
478     }
479   } else {
480     oii  = NULL;
481     ojj  = NULL;
482     oa   = NULL;
483     if (reuse == MAT_INITIAL_MATRIX) {
484       Mat_SeqAIJ* b;
485       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,B);CHKERRQ(ierr);
486       /* hack SeqAIJ */
487       b          = (Mat_SeqAIJ*)((*B)->data);
488       b->free_a  = PETSC_TRUE;
489       b->free_ij = PETSC_TRUE;
490     } else if (reuse == MAT_INPLACE_MATRIX) {
491       Mat T;
492       ierr = MatCreateSeqAIJWithArrays(comm,m,n,dii,djj,da,&T);CHKERRQ(ierr);
493       hypre_CSRMatrixI(hdiag)    = NULL;
494       hypre_CSRMatrixJ(hdiag)    = NULL;
495       hypre_CSRMatrixData(hdiag) = NULL;
496       ierr = MatHeaderReplace(A,&T);CHKERRQ(ierr);
497     }
498   }
499 
500   /* we have to use hypre_Tfree to free the arrays */
501   if (reuse == MAT_INPLACE_MATRIX) {
502     void *ptrs[6] = {dii,djj,da,oii,ojj,oa};
503     const char *names[6] = {"_hypre_csr_dii",
504                             "_hypre_csr_djj",
505                             "_hypre_csr_da",
506                             "_hypre_csr_oii",
507                             "_hypre_csr_ojj",
508                             "_hypre_csr_oa"};
509     for (i=0; i<6; i++) {
510       PetscContainer c;
511 
512       ierr = PetscContainerCreate(comm,&c);CHKERRQ(ierr);
513       ierr = PetscContainerSetPointer(c,ptrs[i]);CHKERRQ(ierr);
514       ierr = PetscContainerSetUserDestroy(c,hypre_array_destroy);CHKERRQ(ierr);
515       ierr = PetscObjectCompose((PetscObject)(*B),names[i],(PetscObject)c);CHKERRQ(ierr);
516       ierr = PetscContainerDestroy(&c);CHKERRQ(ierr);
517     }
518   }
519   PetscFunctionReturn(0);
520 }
521 
522 static PetscErrorCode MatAIJGetParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
523 {
524   hypre_ParCSRMatrix *tA;
525   hypre_CSRMatrix    *hdiag,*hoffd;
526   Mat_SeqAIJ         *diag,*offd;
527   PetscInt           *garray,noffd,dnnz,onnz,*row_starts,*col_starts;
528   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
529   PetscBool          ismpiaij,isseqaij;
530   PetscErrorCode     ierr;
531 
532   PetscFunctionBegin;
533   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
534   ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
535   if (!ismpiaij && !isseqaij) SETERRQ1(comm,PETSC_ERR_SUP,"Unsupported type %s",((PetscObject)A)->type);
536   if (ismpiaij) {
537     Mat_MPIAIJ *a = (Mat_MPIAIJ*)(A->data);
538 
539     diag   = (Mat_SeqAIJ*)a->A->data;
540     offd   = (Mat_SeqAIJ*)a->B->data;
541     garray = a->garray;
542     noffd  = a->B->cmap->N;
543     dnnz   = diag->nz;
544     onnz   = offd->nz;
545   } else {
546     diag   = (Mat_SeqAIJ*)A->data;
547     offd   = NULL;
548     garray = NULL;
549     noffd  = 0;
550     dnnz   = diag->nz;
551     onnz   = 0;
552   }
553 
554   /* create a temporary ParCSR */
555   if (HYPRE_AssumedPartitionCheck()) {
556     PetscMPIInt myid;
557 
558     ierr       = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
559     row_starts = A->rmap->range + myid;
560     col_starts = A->cmap->range + myid;
561   } else {
562     row_starts = A->rmap->range;
563     col_starts = A->cmap->range;
564   }
565   tA = hypre_ParCSRMatrixCreate(comm,A->rmap->N,A->cmap->N,(HYPRE_Int*)row_starts,(HYPRE_Int*)col_starts,noffd,dnnz,onnz);
566   hypre_ParCSRMatrixSetRowStartsOwner(tA,0);
567   hypre_ParCSRMatrixSetColStartsOwner(tA,0);
568 
569   /* set diagonal part */
570   hdiag = hypre_ParCSRMatrixDiag(tA);
571   hypre_CSRMatrixI(hdiag)           = (HYPRE_Int*)diag->i;
572   hypre_CSRMatrixJ(hdiag)           = (HYPRE_Int*)diag->j;
573   hypre_CSRMatrixData(hdiag)        = diag->a;
574   hypre_CSRMatrixNumNonzeros(hdiag) = diag->nz;
575   hypre_CSRMatrixSetRownnz(hdiag);
576   hypre_CSRMatrixSetDataOwner(hdiag,0);
577 
578   /* set offdiagonal part */
579   hoffd = hypre_ParCSRMatrixOffd(tA);
580   if (offd) {
581     hypre_CSRMatrixI(hoffd)           = (HYPRE_Int*)offd->i;
582     hypre_CSRMatrixJ(hoffd)           = (HYPRE_Int*)offd->j;
583     hypre_CSRMatrixData(hoffd)        = offd->a;
584     hypre_CSRMatrixNumNonzeros(hoffd) = offd->nz;
585     hypre_CSRMatrixSetRownnz(hoffd);
586     hypre_CSRMatrixSetDataOwner(hoffd,0);
587     hypre_ParCSRMatrixSetNumNonzeros(tA);
588     hypre_ParCSRMatrixColMapOffd(tA) = (HYPRE_Int*)garray;
589   }
590   *hA = tA;
591   PetscFunctionReturn(0);
592 }
593 
594 static PetscErrorCode MatAIJRestoreParCSR_Private(Mat A, hypre_ParCSRMatrix **hA)
595 {
596   hypre_CSRMatrix    *hdiag,*hoffd;
597 
598   PetscFunctionBegin;
599   hdiag = hypre_ParCSRMatrixDiag(*hA);
600   hoffd = hypre_ParCSRMatrixOffd(*hA);
601   /* set pointers to NULL before destroying tA */
602   hypre_CSRMatrixI(hdiag)           = NULL;
603   hypre_CSRMatrixJ(hdiag)           = NULL;
604   hypre_CSRMatrixData(hdiag)        = NULL;
605   hypre_CSRMatrixI(hoffd)           = NULL;
606   hypre_CSRMatrixJ(hoffd)           = NULL;
607   hypre_CSRMatrixData(hoffd)        = NULL;
608   hypre_ParCSRMatrixColMapOffd(*hA) = NULL;
609   hypre_ParCSRMatrixDestroy(*hA);
610   *hA = NULL;
611   PetscFunctionReturn(0);
612 }
613 
614 /* calls RAP from BoomerAMG:
615    the resulting ParCSR will not own the column and row starts
616    It looks like we don't need to have the diagonal entries
617    ordered first in the rows of the diagonal part
618    for boomerAMGBuildCoarseOperator to work */
619 static PetscErrorCode MatHYPRE_ParCSR_RAP(hypre_ParCSRMatrix *hR, hypre_ParCSRMatrix *hA,hypre_ParCSRMatrix *hP, hypre_ParCSRMatrix **hRAP)
620 {
621   PetscErrorCode ierr;
622   HYPRE_Int      P_owns_col_starts,R_owns_row_starts;
623 
624   PetscFunctionBegin;
625   P_owns_col_starts = hypre_ParCSRMatrixOwnsColStarts(hP);
626   R_owns_row_starts = hypre_ParCSRMatrixOwnsRowStarts(hR);
627   PetscStackCallStandard(hypre_BoomerAMGBuildCoarseOperator,(hR,hA,hP,hRAP));
628   PetscStackCallStandard(hypre_ParCSRMatrixSetNumNonzeros,(*hRAP));
629   /* hypre_BoomerAMGBuildCoarseOperator steals the col_starts from P and the row_starts from R */
630   hypre_ParCSRMatrixSetRowStartsOwner(*hRAP,0);
631   hypre_ParCSRMatrixSetColStartsOwner(*hRAP,0);
632   if (P_owns_col_starts) hypre_ParCSRMatrixSetColStartsOwner(hP,1);
633   if (R_owns_row_starts) hypre_ParCSRMatrixSetRowStartsOwner(hR,1);
634   PetscFunctionReturn(0);
635 }
636 
637 static PetscErrorCode MatPtAPNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat P,Mat C)
638 {
639   Mat                B;
640   hypre_ParCSRMatrix *hA,*hP,*hPtAP;
641   PetscErrorCode     ierr;
642 
643   PetscFunctionBegin;
644   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
645   ierr = MatAIJGetParCSR_Private(P,&hP);CHKERRQ(ierr);
646   ierr = MatHYPRE_ParCSR_RAP(hP,hA,hP,&hPtAP);CHKERRQ(ierr);
647   ierr = MatCreateFromParCSR(hPtAP,MATAIJ,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
648   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
649   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
650   ierr = MatAIJRestoreParCSR_Private(P,&hP);CHKERRQ(ierr);
651   PetscFunctionReturn(0);
652 }
653 
654 PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat P,PetscReal fill,Mat *C)
655 {
656   PetscErrorCode ierr;
657 
658   PetscFunctionBegin;
659   ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
660   ierr                   = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
661   (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_AIJ_wHYPRE;
662   PetscFunctionReturn(0);
663 }
664 
665 static PetscErrorCode MatPtAPNumeric_AIJ_HYPRE(Mat A,Mat P,Mat C)
666 {
667   Mat                B;
668   Mat_HYPRE          *hP;
669   hypre_ParCSRMatrix *hA = NULL,*Pparcsr,*ptapparcsr;
670   HYPRE_Int          type;
671   MPI_Comm           comm = PetscObjectComm((PetscObject)A);
672   PetscBool          ishypre;
673   PetscErrorCode     ierr;
674 
675   PetscFunctionBegin;
676   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
677   if (!ishypre) SETERRQ1(comm,PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
678   hP = (Mat_HYPRE*)P->data;
679   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
680   if (type != HYPRE_PARCSR) SETERRQ(comm,PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
681   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
682 
683   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
684   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,hA,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
685   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
686 
687   /* create temporary matrix and merge to C */
688   ierr = MatCreateFromParCSR(ptapparcsr,((PetscObject)C)->type_name,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
689   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
690   PetscFunctionReturn(0);
691 }
692 
693 static PetscErrorCode MatPtAP_AIJ_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
694 {
695   PetscErrorCode ierr;
696 
697   PetscFunctionBegin;
698   if (scall == MAT_INITIAL_MATRIX) {
699     const char *deft = MATAIJ;
700     char       type[256];
701     PetscBool  flg;
702 
703     ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr);
704     ierr = PetscOptionsFList("-matptap_hypre_outtype","Matrix type",NULL,MatList,deft,type,256,&flg);CHKERRQ(ierr);
705     ierr = PetscOptionsEnd();CHKERRQ(ierr);
706     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
707     ierr = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
708     if (flg) {
709       ierr = MatSetType(*C,type);CHKERRQ(ierr);
710     } else {
711       ierr = MatSetType(*C,deft);CHKERRQ(ierr);
712     }
713     (*C)->ops->ptapnumeric = MatPtAPNumeric_AIJ_HYPRE;
714     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
715   }
716   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
717   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
718   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
719   PetscFunctionReturn(0);
720 }
721 
722 static PetscErrorCode MatPtAPNumeric_HYPRE_HYPRE(Mat A,Mat P,Mat C)
723 {
724   Mat                B;
725   hypre_ParCSRMatrix *Aparcsr,*Pparcsr,*ptapparcsr;
726   Mat_HYPRE          *hA,*hP;
727   PetscBool          ishypre;
728   HYPRE_Int          type;
729   PetscErrorCode     ierr;
730 
731   PetscFunctionBegin;
732   ierr = PetscObjectTypeCompare((PetscObject)P,MATHYPRE,&ishypre);CHKERRQ(ierr);
733   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)P),PETSC_ERR_USER,"P should be of type %s",MATHYPRE);
734   ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr);
735   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE);
736   hA = (Mat_HYPRE*)A->data;
737   hP = (Mat_HYPRE*)P->data;
738   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
739   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
740   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hP->ij,&type));
741   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)P),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
742   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
743   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hP->ij,(void**)&Pparcsr));
744   ierr = MatHYPRE_ParCSR_RAP(Pparcsr,Aparcsr,Pparcsr,&ptapparcsr);CHKERRQ(ierr);
745   ierr = MatCreateFromParCSR(ptapparcsr,MATHYPRE,PETSC_OWN_POINTER,&B);CHKERRQ(ierr);
746   ierr = MatHeaderMerge(C,&B);CHKERRQ(ierr);
747   PetscFunctionReturn(0);
748 }
749 
750 static PetscErrorCode MatPtAP_HYPRE_HYPRE(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
751 {
752   PetscErrorCode ierr;
753 
754   PetscFunctionBegin;
755   if (scall == MAT_INITIAL_MATRIX) {
756     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
757     ierr                   = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
758     ierr                   = MatSetType(*C,MATHYPRE);CHKERRQ(ierr);
759     (*C)->ops->ptapnumeric = MatPtAPNumeric_HYPRE_HYPRE;
760     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
761   }
762   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
763   ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr);
764   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
765   PetscFunctionReturn(0);
766 }
767 
768 /* calls hypre_ParMatmul
769    hypre_ParMatMul uses hypre_ParMatrixCreate with the communicator of hA
770    hypre_ParMatrixCreate does not duplicate the communicator
771    It looks like we don't need to have the diagonal entries
772    ordered first in the rows of the diagonal part
773    for boomerAMGBuildCoarseOperator to work */
774 static PetscErrorCode MatHYPRE_ParCSR_MatMatMult(hypre_ParCSRMatrix *hA, hypre_ParCSRMatrix *hB, hypre_ParCSRMatrix **hAB)
775 {
776   PetscFunctionBegin;
777   PetscStackPush("hypre_ParMatmul");
778   *hAB = hypre_ParMatmul(hA,hB);
779   PetscStackPop;
780   PetscFunctionReturn(0);
781 }
782 
783 static PetscErrorCode MatMatMultNumeric_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C)
784 {
785   Mat                D;
786   hypre_ParCSRMatrix *hA,*hB,*hAB = NULL;
787   PetscErrorCode     ierr;
788 
789   PetscFunctionBegin;
790   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
791   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
792   ierr = MatHYPRE_ParCSR_MatMatMult(hA,hB,&hAB);CHKERRQ(ierr);
793   ierr = MatCreateFromParCSR(hAB,MATAIJ,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
794   ierr = MatHeaderMerge(C,&D);CHKERRQ(ierr);
795   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
796   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
797   PetscFunctionReturn(0);
798 }
799 
800 PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat A,Mat B,PetscReal fill,Mat *C)
801 {
802   PetscErrorCode ierr;
803 
804   PetscFunctionBegin;
805   ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
806   ierr                      = MatSetType(*C,MATAIJ);CHKERRQ(ierr);
807   (*C)->ops->matmultnumeric = MatMatMultNumeric_AIJ_AIJ_wHYPRE;
808   PetscFunctionReturn(0);
809 }
810 
811 static PetscErrorCode MatMatMultNumeric_HYPRE_HYPRE(Mat A,Mat B,Mat C)
812 {
813   Mat                D;
814   hypre_ParCSRMatrix *Aparcsr,*Bparcsr,*ABparcsr = NULL;
815   Mat_HYPRE          *hA,*hB;
816   PetscBool          ishypre;
817   HYPRE_Int          type;
818   PetscErrorCode     ierr;
819 
820   PetscFunctionBegin;
821   ierr = PetscObjectTypeCompare((PetscObject)B,MATHYPRE,&ishypre);CHKERRQ(ierr);
822   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)B),PETSC_ERR_USER,"B should be of type %s",MATHYPRE);
823   ierr = PetscObjectTypeCompare((PetscObject)A,MATHYPRE,&ishypre);CHKERRQ(ierr);
824   if (!ishypre) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"A should be of type %s",MATHYPRE);
825   hA = (Mat_HYPRE*)A->data;
826   hB = (Mat_HYPRE*)B->data;
827   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
828   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
829   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hB->ij,&type));
830   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_SUP,"Only HYPRE_PARCSR is supported");
831   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&Aparcsr));
832   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hB->ij,(void**)&Bparcsr));
833   ierr = MatHYPRE_ParCSR_MatMatMult(Aparcsr,Bparcsr,&ABparcsr);CHKERRQ(ierr);
834   ierr = MatCreateFromParCSR(ABparcsr,MATHYPRE,PETSC_OWN_POINTER,&D);CHKERRQ(ierr);
835   /* need to use HeaderReplace because HeaderMerge messes up with the communicator */
836   ierr = MatHeaderReplace(C,&D);CHKERRQ(ierr);
837   C->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE;
838   PetscFunctionReturn(0);
839 }
840 
841 static PetscErrorCode MatMatMult_HYPRE_HYPRE(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
842 {
843   PetscErrorCode ierr;
844 
845   PetscFunctionBegin;
846   if (scall == MAT_INITIAL_MATRIX) {
847     ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
848     ierr                      = MatCreate(PetscObjectComm((PetscObject)A),C);CHKERRQ(ierr);
849     ierr                      = MatSetType(*C,MATHYPRE);CHKERRQ(ierr);
850     (*C)->ops->matmultnumeric = MatMatMultNumeric_HYPRE_HYPRE;
851     ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr);
852   }
853   ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
854   ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr);
855   ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr);
856   PetscFunctionReturn(0);
857 }
858 
859 PETSC_INTERN PetscErrorCode MatTransposeMatMatMultNumeric_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,Mat D)
860 {
861   Mat                E;
862   hypre_ParCSRMatrix *hA,*hB,*hC,*hABC;
863   PetscErrorCode     ierr;
864 
865   PetscFunctionBegin;
866   ierr = MatAIJGetParCSR_Private(A,&hA);CHKERRQ(ierr);
867   ierr = MatAIJGetParCSR_Private(B,&hB);CHKERRQ(ierr);
868   ierr = MatAIJGetParCSR_Private(C,&hC);CHKERRQ(ierr);
869   ierr = MatHYPRE_ParCSR_RAP(hA,hB,hC,&hABC);CHKERRQ(ierr);
870   ierr = MatCreateFromParCSR(hABC,MATAIJ,PETSC_OWN_POINTER,&E);CHKERRQ(ierr);
871   ierr = MatHeaderMerge(D,&E);CHKERRQ(ierr);
872   ierr = MatAIJRestoreParCSR_Private(A,&hA);CHKERRQ(ierr);
873   ierr = MatAIJRestoreParCSR_Private(B,&hB);CHKERRQ(ierr);
874   ierr = MatAIJRestoreParCSR_Private(C,&hC);CHKERRQ(ierr);
875   PetscFunctionReturn(0);
876 }
877 
878 PETSC_INTERN PetscErrorCode MatTransposeMatMatMultSymbolic_AIJ_AIJ_AIJ_wHYPRE(Mat A,Mat B,Mat C,PetscReal fill,Mat *D)
879 {
880   PetscErrorCode ierr;
881 
882   PetscFunctionBegin;
883   ierr = MatCreate(PetscObjectComm((PetscObject)A),D);CHKERRQ(ierr);
884   ierr = MatSetType(*D,MATAIJ);CHKERRQ(ierr);
885   PetscFunctionReturn(0);
886 }
887 
888 static PetscErrorCode MatMultTranspose_HYPRE(Mat A, Vec x, Vec y)
889 {
890   PetscErrorCode ierr;
891 
892   PetscFunctionBegin;
893   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_TRUE);CHKERRQ(ierr);
894   PetscFunctionReturn(0);
895 }
896 
897 static PetscErrorCode MatMult_HYPRE(Mat A, Vec x, Vec y)
898 {
899   PetscErrorCode ierr;
900 
901   PetscFunctionBegin;
902   ierr = MatHYPRE_MultKernel_Private(A,x,y,PETSC_FALSE);CHKERRQ(ierr);
903   PetscFunctionReturn(0);
904 }
905 
906 static PetscErrorCode MatHYPRE_MultKernel_Private(Mat A, Vec x, Vec y, PetscBool trans)
907 {
908   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
909   hypre_ParCSRMatrix *parcsr;
910   hypre_ParVector    *hx,*hy;
911   PetscScalar        *ax,*ay,*sax,*say;
912   PetscErrorCode     ierr;
913 
914   PetscFunctionBegin;
915   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)&parcsr));
916   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->x,(void**)&hx));
917   PetscStackCallStandard(HYPRE_IJVectorGetObject,(hA->b,(void**)&hy));
918   ierr = VecGetArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
919   ierr = VecGetArray(y,&ay);CHKERRQ(ierr);
920   if (trans) {
921     VecHYPRE_ParVectorReplacePointer(hA->x,ay,say);
922     VecHYPRE_ParVectorReplacePointer(hA->b,ax,sax);
923     hypre_ParCSRMatrixMatvecT(1.,parcsr,hy,0.,hx);
924     VecHYPRE_ParVectorReplacePointer(hA->x,say,ay);
925     VecHYPRE_ParVectorReplacePointer(hA->b,sax,ax);
926   } else {
927     VecHYPRE_ParVectorReplacePointer(hA->x,ax,sax);
928     VecHYPRE_ParVectorReplacePointer(hA->b,ay,say);
929     hypre_ParCSRMatrixMatvec(1.,parcsr,hx,0.,hy);
930     VecHYPRE_ParVectorReplacePointer(hA->x,sax,ax);
931     VecHYPRE_ParVectorReplacePointer(hA->b,say,ay);
932   }
933   ierr = VecRestoreArrayRead(x,(const PetscScalar**)&ax);CHKERRQ(ierr);
934   ierr = VecRestoreArray(y,&ay);CHKERRQ(ierr);
935   PetscFunctionReturn(0);
936 }
937 
938 static PetscErrorCode MatDestroy_HYPRE(Mat A)
939 {
940   Mat_HYPRE      *hA = (Mat_HYPRE*)A->data;
941   PetscErrorCode ierr;
942 
943   PetscFunctionBegin;
944   if (hA->x) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->x));
945   if (hA->b) PetscStackCallStandard(HYPRE_IJVectorDestroy,(hA->b));
946   if (hA->ij) {
947     if (!hA->inner_free) hypre_IJMatrixObject(hA->ij) = NULL;
948     PetscStackCallStandard(HYPRE_IJMatrixDestroy,(hA->ij));
949   }
950   if (hA->comm) { ierr = MPI_Comm_free(&hA->comm);CHKERRQ(ierr); }
951 
952   ierr = MatStashDestroy_Private(&A->stash);CHKERRQ(ierr);
953 
954   ierr = PetscFree(hA->array);CHKERRQ(ierr);
955 
956   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_aij_C",NULL);CHKERRQ(ierr);
957   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_hypre_is_C",NULL);CHKERRQ(ierr);
958   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_seqaij_hypre_C",NULL);CHKERRQ(ierr);
959   ierr = PetscObjectComposeFunction((PetscObject)A,"MatPtAP_mpiaij_hypre_C",NULL);CHKERRQ(ierr);
960   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPRESetPreallocation_C",NULL);CHKERRQ(ierr);
961   ierr = PetscObjectComposeFunction((PetscObject)A,"MatHYPREGetParCSR_C",NULL);CHKERRQ(ierr);
962   ierr = PetscFree(A->data);CHKERRQ(ierr);
963   PetscFunctionReturn(0);
964 }
965 
966 static PetscErrorCode MatSetUp_HYPRE(Mat A)
967 {
968   PetscErrorCode ierr;
969 
970   PetscFunctionBegin;
971   ierr = MatHYPRESetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr);
972   PetscFunctionReturn(0);
973 }
974 
975 
976 static PetscErrorCode MatAssemblyEnd_HYPRE(Mat A, MatAssemblyType mode)
977 {
978   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
979   Vec                x,b;
980   PetscMPIInt        n;
981   PetscInt           i,j,rstart,ncols,flg;
982   PetscInt           *row,*col;
983   PetscScalar        *val;
984   PetscErrorCode     ierr;
985 
986   PetscFunctionBegin;
987   if (mode == MAT_FLUSH_ASSEMBLY) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_FLUSH_ASSEMBLY currently not supported with MATHYPRE");
988 
989   if (!A->nooffprocentries) {
990     while (1) {
991       ierr = MatStashScatterGetMesg_Private(&A->stash,&n,&row,&col,&val,&flg);CHKERRQ(ierr);
992       if (!flg) break;
993 
994       for (i=0; i<n; ) {
995         /* Now identify the consecutive vals belonging to the same row */
996         for (j=i,rstart=row[j]; j<n; j++) {
997           if (row[j] != rstart) break;
998         }
999         if (j < n) ncols = j-i;
1000         else       ncols = n-i;
1001         /* Now assemble all these values with a single function call */
1002         ierr = MatSetValues_HYPRE(A,1,row+i,ncols,col+i,val+i,A->insertmode);CHKERRQ(ierr);
1003 
1004         i = j;
1005       }
1006     }
1007     ierr = MatStashScatterEnd_Private(&A->stash);CHKERRQ(ierr);
1008   }
1009 
1010   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(hA->ij));
1011   if (hA->x) PetscFunctionReturn(0);
1012   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1013   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1014   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->cmap->n,A->cmap->N,NULL,&x);CHKERRQ(ierr);
1015   ierr = VecCreateMPIWithArray(PetscObjectComm((PetscObject)A),1,A->rmap->n,A->rmap->N,NULL,&b);CHKERRQ(ierr);
1016   ierr = VecHYPRE_IJVectorCreate(x,&hA->x);CHKERRQ(ierr);
1017   ierr = VecHYPRE_IJVectorCreate(b,&hA->b);CHKERRQ(ierr);
1018   ierr = VecDestroy(&x);CHKERRQ(ierr);
1019   ierr = VecDestroy(&b);CHKERRQ(ierr);
1020   PetscFunctionReturn(0);
1021 }
1022 
1023 static PetscErrorCode MatGetArray_HYPRE(Mat A, PetscInt size, void **array)
1024 {
1025   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
1026   PetscErrorCode     ierr;
1027 
1028   PetscFunctionBegin;
1029   if (!hA->available) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Temporary space is in use");
1030 
1031   if (hA->size >= size) *array = hA->array;
1032   else {
1033     ierr = PetscFree(hA->array);CHKERRQ(ierr);
1034     hA->size = size;
1035     ierr = PetscMalloc(hA->size,&hA->array);CHKERRQ(ierr);
1036     *array = hA->array;
1037   }
1038 
1039   hA->available = PETSC_FALSE;
1040   PetscFunctionReturn(0);
1041 }
1042 
1043 static PetscErrorCode MatRestoreArray_HYPRE(Mat A, PetscInt size, void **array)
1044 {
1045   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
1046 
1047   PetscFunctionBegin;
1048   *array = NULL;
1049   size = 0;
1050   hA->available = PETSC_TRUE;
1051   PetscFunctionReturn(0);
1052 }
1053 
1054 
1055 PetscErrorCode MatSetValues_HYPRE(Mat A, PetscInt nr, const PetscInt rows[], PetscInt nc, const PetscInt cols[], const PetscScalar v[], InsertMode ins)
1056 {
1057   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
1058   PetscScalar        *vals = (PetscScalar *)v;
1059   PetscScalar        *sscr;
1060   PetscInt           *cscr[2];
1061   PetscInt           i,nzc;
1062   void               *array;
1063   PetscErrorCode     ierr;
1064 
1065   PetscFunctionBegin;
1066   ierr = MatGetArray_HYPRE(A,sizeof(PetscInt)*(2*nc)+sizeof(PetscScalar)*nc*nr,&array);CHKERRQ(ierr);
1067   cscr[0] = (PetscInt*)array;
1068   cscr[1] = ((PetscInt*)array)+nc;
1069   sscr = (PetscScalar*)(((PetscInt*)array)+nc*2);
1070   for (i=0,nzc=0;i<nc;i++) {
1071     if (cols[i] >= 0) {
1072       cscr[0][nzc  ] = cols[i];
1073       cscr[1][nzc++] = i;
1074     }
1075   }
1076   if (!nzc) {
1077     ierr = MatRestoreArray_HYPRE(A,sizeof(PetscInt)*(2*nc)+sizeof(PetscScalar)*nc*nr,&array);CHKERRQ(ierr);
1078     PetscFunctionReturn(0);
1079   }
1080 
1081   if (ins == ADD_VALUES) {
1082     for (i=0;i<nr;i++) {
1083       if (rows[i] >= 0 && nzc) {
1084         PetscInt j;
1085         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
1086         PetscStackCallStandard(HYPRE_IJMatrixAddToValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
1087       }
1088       vals += nc;
1089     }
1090   } else { /* INSERT_VALUES */
1091 
1092     PetscInt rst,ren;
1093     ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
1094 
1095     for (i=0;i<nr;i++) {
1096       if (rows[i] >= 0 && nzc) {
1097         PetscInt j;
1098         for (j=0;j<nzc;j++) sscr[j] = vals[cscr[1][j]];
1099         /* nonlocal values */
1100         if (rows[i] < rst || rows[i] >= ren) { ierr = MatStashValuesRow_Private(&A->stash,rows[i],nzc,cscr[0],sscr, PETSC_FALSE);CHKERRQ(ierr); }
1101         /* local values */
1102         else PetscStackCallStandard(HYPRE_IJMatrixSetValues,(hA->ij,1,&nzc,(HYPRE_Int*)(rows+i),cscr[0],sscr));
1103       }
1104       vals += nc;
1105     }
1106   }
1107 
1108   ierr = MatRestoreArray_HYPRE(A,sizeof(PetscInt)*(2*nc)+sizeof(PetscScalar)*nc*nr,&array);CHKERRQ(ierr);
1109   PetscFunctionReturn(0);
1110 }
1111 
1112 static PetscErrorCode MatHYPRESetPreallocation_HYPRE(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1113 {
1114   Mat_HYPRE          *hA = (Mat_HYPRE*)A->data;
1115   HYPRE_Int          *hdnnz,*honnz;
1116   PetscInt           i,rs,re,cs,ce,bs;
1117   PetscMPIInt        size;
1118   PetscErrorCode     ierr;
1119 
1120   PetscFunctionBegin;
1121   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
1122   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1123   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1124   rs   = A->rmap->rstart;
1125   re   = A->rmap->rend;
1126   cs   = A->cmap->rstart;
1127   ce   = A->cmap->rend;
1128   if (!hA->ij) {
1129     PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rs,re-1,cs,ce-1,&hA->ij));
1130     PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1131   } else {
1132     HYPRE_Int hrs,hre,hcs,hce;
1133     PetscStackCallStandard(HYPRE_IJMatrixGetLocalRange,(hA->ij,&hrs,&hre,&hcs,&hce));
1134     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);
1135     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);
1136   }
1137   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1138 
1139   if (dnz == PETSC_DEFAULT || dnz == PETSC_DECIDE) dnz = 10*bs;
1140   if (onz == PETSC_DEFAULT || onz == PETSC_DECIDE) onz = 10*bs;
1141 
1142   if (!dnnz) {
1143     ierr = PetscMalloc1(A->rmap->n,&hdnnz);CHKERRQ(ierr);
1144     for (i=0;i<A->rmap->n;i++) hdnnz[i] = dnz;
1145   } else {
1146     hdnnz = (HYPRE_Int*)dnnz;
1147   }
1148   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1149   if (size > 1) {
1150     if (!onnz) {
1151       ierr = PetscMalloc1(A->rmap->n,&honnz);CHKERRQ(ierr);
1152       for (i=0;i<A->rmap->n;i++) honnz[i] = onz;
1153     } else {
1154       honnz = (HYPRE_Int*)onnz;
1155     }
1156     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(hA->ij,hdnnz,honnz));
1157   } else {
1158     honnz = NULL;
1159     PetscStackCallStandard(HYPRE_IJMatrixSetRowSizes,(hA->ij,hdnnz));
1160   }
1161   if (!dnnz) {
1162     ierr = PetscFree(hdnnz);CHKERRQ(ierr);
1163   }
1164   if (!onnz && honnz) {
1165     ierr = PetscFree(honnz);CHKERRQ(ierr);
1166   }
1167   A->preallocated = PETSC_TRUE;
1168 
1169   /* SetDiagOffdSizes sets hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0 */
1170   {
1171     hypre_AuxParCSRMatrix *aux_matrix;
1172     aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(hA->ij);
1173     hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 1;
1174   }
1175   PetscFunctionReturn(0);
1176 }
1177 
1178 /*@C
1179    MatHYPRESetPreallocation - Preallocates memory for a sparse parallel matrix in HYPRE IJ format
1180 
1181    Collective on Mat
1182 
1183    Input Parameters:
1184 +  A - the matrix
1185 .  dnz  - number of nonzeros per row in DIAGONAL portion of local submatrix
1186           (same value is used for all local rows)
1187 .  dnnz - array containing the number of nonzeros in the various rows of the
1188           DIAGONAL portion of the local submatrix (possibly different for each row)
1189           or NULL (PETSC_NULL_INTEGER in Fortran), if d_nz is used to specify the nonzero structure.
1190           The size of this array is equal to the number of local rows, i.e 'm'.
1191           For matrices that will be factored, you must leave room for (and set)
1192           the diagonal entry even if it is zero.
1193 .  onz  - number of nonzeros per row in the OFF-DIAGONAL portion of local
1194           submatrix (same value is used for all local rows).
1195 -  onnz - array containing the number of nonzeros in the various rows of the
1196           OFF-DIAGONAL portion of the local submatrix (possibly different for
1197           each row) or NULL (PETSC_NULL_INTEGER in Fortran), if o_nz is used to specify the nonzero
1198           structure. The size of this array is equal to the number
1199           of local rows, i.e 'm'.
1200 
1201    Notes: If the *nnz parameter is given then the *nz parameter is ignored; for sequential matrices, onz and onnz are ignored.
1202 
1203    Level: intermediate
1204 
1205 .keywords: matrix, aij, compressed row, sparse, parallel
1206 
1207 .seealso: MatCreate(), MatMPIAIJSetPreallocation, MATHYPRE
1208 @*/
1209 PetscErrorCode MatHYPRESetPreallocation(Mat A, PetscInt dnz, const PetscInt dnnz[], PetscInt onz, const PetscInt onnz[])
1210 {
1211   PetscErrorCode ierr;
1212 
1213   PetscFunctionBegin;
1214   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1215   PetscValidType(A,1);
1216   ierr = PetscTryMethod(A,"MatHYPRESetPreallocation_C",(Mat,PetscInt,const PetscInt[],PetscInt,const PetscInt[]),(A,dnz,dnnz,onz,onnz));CHKERRQ(ierr);
1217   PetscFunctionReturn(0);
1218 }
1219 
1220 /*
1221    MatCreateFromParCSR - Creates a matrix from a hypre_ParCSRMatrix
1222 
1223    Collective
1224 
1225    Input Parameters:
1226 +  vparcsr  - the pointer to the hypre_ParCSRMatrix
1227 .  mtype    - matrix type to be created. Currently MATAIJ, MATIS and MATHYPRE are supported.
1228 -  copymode - PETSc copying options
1229 
1230    Output Parameter:
1231 .  A  - the matrix
1232 
1233    Level: intermediate
1234 
1235 .seealso: MatHYPRE, PetscCopyMode
1236 */
1237 PETSC_EXTERN PetscErrorCode MatCreateFromParCSR(hypre_ParCSRMatrix *vparcsr, MatType mtype, PetscCopyMode copymode, Mat* A)
1238 {
1239   Mat                   T;
1240   Mat_HYPRE             *hA;
1241   hypre_ParCSRMatrix    *parcsr;
1242   MPI_Comm              comm;
1243   PetscInt              rstart,rend,cstart,cend,M,N;
1244   PetscBool             isseqaij,ismpiaij,isaij,ishyp,isis;
1245   PetscErrorCode        ierr;
1246 
1247   PetscFunctionBegin;
1248   parcsr = (hypre_ParCSRMatrix *)vparcsr;
1249   comm   = hypre_ParCSRMatrixComm(parcsr);
1250   ierr   = PetscStrcmp(mtype,MATSEQAIJ,&isseqaij);CHKERRQ(ierr);
1251   ierr   = PetscStrcmp(mtype,MATMPIAIJ,&ismpiaij);CHKERRQ(ierr);
1252   ierr   = PetscStrcmp(mtype,MATAIJ,&isaij);CHKERRQ(ierr);
1253   ierr   = PetscStrcmp(mtype,MATHYPRE,&ishyp);CHKERRQ(ierr);
1254   ierr   = PetscStrcmp(mtype,MATIS,&isis);CHKERRQ(ierr);
1255   isaij  = (PetscBool)(isseqaij || ismpiaij || isaij);
1256   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);
1257   if (ishyp && copymode == PETSC_COPY_VALUES) SETERRQ(comm,PETSC_ERR_SUP,"Unsupported copymode PETSC_COPY_VALUES");
1258 
1259   /* access ParCSRMatrix */
1260   rstart = hypre_ParCSRMatrixFirstRowIndex(parcsr);
1261   rend   = hypre_ParCSRMatrixLastRowIndex(parcsr);
1262   cstart = hypre_ParCSRMatrixFirstColDiag(parcsr);
1263   cend   = hypre_ParCSRMatrixLastColDiag(parcsr);
1264   M      = hypre_ParCSRMatrixGlobalNumRows(parcsr);
1265   N      = hypre_ParCSRMatrixGlobalNumCols(parcsr);
1266 
1267   /* fix for empty local rows/columns */
1268   if (rend < rstart) rend = rstart;
1269   if (cend < cstart) cend = cstart;
1270 
1271   /* create PETSc matrix with MatHYPRE */
1272   ierr = MatCreate(comm,&T);CHKERRQ(ierr);
1273   ierr = MatSetSizes(T,rend-rstart+1,cend-cstart+1,M,N);CHKERRQ(ierr);
1274   ierr = MatSetType(T,MATHYPRE);CHKERRQ(ierr);
1275   hA   = (Mat_HYPRE*)(T->data);
1276 
1277   /* create HYPRE_IJMatrix */
1278   PetscStackCallStandard(HYPRE_IJMatrixCreate,(hA->comm,rstart,rend-1,cstart,cend-1,&hA->ij));
1279 
1280   /* set ParCSR object */
1281   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(hA->ij,HYPRE_PARCSR));
1282   hypre_IJMatrixObject(hA->ij) = parcsr;
1283   T->preallocated = PETSC_TRUE;
1284 
1285   /* set assembled flag */
1286   hypre_IJMatrixAssembleFlag(hA->ij) = 1;
1287   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(hA->ij));
1288   if (ishyp) {
1289     PetscMPIInt myid = 0;
1290 
1291     /* make sure we always have row_starts and col_starts available */
1292     if (HYPRE_AssumedPartitionCheck()) {
1293       ierr = MPI_Comm_rank(comm,&myid);CHKERRQ(ierr);
1294     }
1295     if (!hypre_ParCSRMatrixOwnsColStarts(parcsr)) {
1296       PetscLayout map;
1297 
1298       ierr = MatGetLayouts(T,NULL,&map);CHKERRQ(ierr);
1299       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1300       hypre_ParCSRMatrixColStarts(parcsr) = (HYPRE_Int*)(map->range + myid);
1301     }
1302     if (!hypre_ParCSRMatrixOwnsRowStarts(parcsr)) {
1303       PetscLayout map;
1304 
1305       ierr = MatGetLayouts(T,&map,NULL);CHKERRQ(ierr);
1306       ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1307       hypre_ParCSRMatrixRowStarts(parcsr) = (HYPRE_Int*)(map->range + myid);
1308     }
1309     /* prevent from freeing the pointer */
1310     if (copymode == PETSC_USE_POINTER) hA->inner_free = PETSC_FALSE;
1311     *A   = T;
1312     ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1313     ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1314   } else if (isaij) {
1315     if (copymode != PETSC_OWN_POINTER) {
1316       /* prevent from freeing the pointer */
1317       hA->inner_free = PETSC_FALSE;
1318       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
1319       ierr = MatDestroy(&T);CHKERRQ(ierr);
1320     } else { /* AIJ return type with PETSC_OWN_POINTER */
1321       ierr = MatConvert_HYPRE_AIJ(T,MATAIJ,MAT_INPLACE_MATRIX,&T);CHKERRQ(ierr);
1322       *A   = T;
1323     }
1324   } else if (isis) {
1325     ierr = MatConvert_HYPRE_IS(T,MATIS,MAT_INITIAL_MATRIX,A);CHKERRQ(ierr);
1326     if (copymode != PETSC_OWN_POINTER) hA->inner_free = PETSC_FALSE;
1327     ierr = MatDestroy(&T);CHKERRQ(ierr);
1328   }
1329   PetscFunctionReturn(0);
1330 }
1331 
1332 static PetscErrorCode MatHYPREGetParCSR_HYPRE(Mat A, hypre_ParCSRMatrix **parcsr)
1333 {
1334   Mat_HYPRE*            hA = (Mat_HYPRE*)A->data;
1335   HYPRE_Int             type;
1336   PetscErrorCode        ierr;
1337 
1338   PetscFunctionBegin;
1339   if (!hA->ij) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"HYPRE_IJMatrix not present");
1340   PetscStackCallStandard(HYPRE_IJMatrixGetObjectType,(hA->ij,&type));
1341   if (type != HYPRE_PARCSR) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"HYPRE_IJMatrix is not of type HYPRE_PARCSR");
1342   PetscStackCallStandard(HYPRE_IJMatrixGetObject,(hA->ij,(void**)parcsr));
1343   PetscFunctionReturn(0);
1344 }
1345 
1346 /*
1347    MatHYPREGetParCSR - Gets the pointer to the ParCSR matrix
1348 
1349    Not collective
1350 
1351    Input Parameters:
1352 +  A  - the MATHYPRE object
1353 
1354    Output Parameter:
1355 .  parcsr  - the pointer to the hypre_ParCSRMatrix
1356 
1357    Level: intermediate
1358 
1359 .seealso: MatHYPRE, PetscCopyMode
1360 */
1361 PetscErrorCode MatHYPREGetParCSR(Mat A, hypre_ParCSRMatrix **parcsr)
1362 {
1363   PetscErrorCode ierr;
1364 
1365   PetscFunctionBegin;
1366   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1367   PetscValidType(A,1);
1368   ierr = PetscUseMethod(A,"MatHYPREGetParCSR_C",(Mat,hypre_ParCSRMatrix**),(A,parcsr));CHKERRQ(ierr);
1369   PetscFunctionReturn(0);
1370 }
1371 
1372 static PetscErrorCode MatMissingDiagonal_HYPRE(Mat A, PetscBool *missing, PetscInt *dd)
1373 {
1374   hypre_ParCSRMatrix *parcsr;
1375   hypre_CSRMatrix    *ha;
1376   PetscInt           rst;
1377   PetscErrorCode     ierr;
1378 
1379   PetscFunctionBegin;
1380   if (A->rmap->n != A->cmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented with non-square diagonal blocks");
1381   ierr = MatGetOwnershipRange(A,&rst,NULL);CHKERRQ(ierr);
1382   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1383   if (missing) *missing = PETSC_FALSE;
1384   if (dd) *dd = -1;
1385   ha = hypre_ParCSRMatrixDiag(parcsr);
1386   if (ha) {
1387     PetscInt  size,i;
1388     HYPRE_Int *ii,*jj;
1389 
1390     size = hypre_CSRMatrixNumRows(ha);
1391     ii   = hypre_CSRMatrixI(ha);
1392     jj   = hypre_CSRMatrixJ(ha);
1393     for (i = 0; i < size; i++) {
1394       PetscInt  j;
1395       PetscBool found = PETSC_FALSE;
1396 
1397       for (j = ii[i]; j < ii[i+1] && !found; j++)
1398         found = (jj[j] == i) ? PETSC_TRUE : PETSC_FALSE;
1399 
1400       if (!found) {
1401         PetscInfo1(A,"Matrix is missing local diagonal entry %D\n",i);
1402         if (missing) *missing = PETSC_TRUE;
1403         if (dd) *dd = i+rst;
1404         PetscFunctionReturn(0);
1405       }
1406     }
1407     if (!size) {
1408       PetscInfo(A,"Matrix has no diagonal entries therefore is missing diagonal\n");
1409       if (missing) *missing = PETSC_TRUE;
1410       if (dd) *dd = rst;
1411     }
1412   } else {
1413     PetscInfo(A,"Matrix has no diagonal entries therefore is missing diagonal\n");
1414     if (missing) *missing = PETSC_TRUE;
1415     if (dd) *dd = rst;
1416   }
1417   PetscFunctionReturn(0);
1418 }
1419 
1420 static PetscErrorCode MatScale_HYPRE(Mat A, PetscScalar s)
1421 {
1422   hypre_ParCSRMatrix *parcsr;
1423   hypre_CSRMatrix    *ha;
1424   PetscErrorCode     ierr;
1425 
1426   PetscFunctionBegin;
1427   ierr  = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1428   /* diagonal part */
1429   ha = hypre_ParCSRMatrixDiag(parcsr);
1430   if (ha) {
1431     PetscInt    size,i;
1432     HYPRE_Int   *ii;
1433     PetscScalar *a;
1434 
1435     size = hypre_CSRMatrixNumRows(ha);
1436     a    = hypre_CSRMatrixData(ha);
1437     ii   = hypre_CSRMatrixI(ha);
1438     for (i = 0; i < ii[size]; i++) a[i] *= s;
1439   }
1440   /* offdiagonal part */
1441   ha = hypre_ParCSRMatrixOffd(parcsr);
1442   if (ha) {
1443     PetscInt    size,i;
1444     HYPRE_Int   *ii;
1445     PetscScalar *a;
1446 
1447     size = hypre_CSRMatrixNumRows(ha);
1448     a    = hypre_CSRMatrixData(ha);
1449     ii   = hypre_CSRMatrixI(ha);
1450     for (i = 0; i < ii[size]; i++) a[i] *= s;
1451   }
1452   PetscFunctionReturn(0);
1453 }
1454 
1455 static PetscErrorCode MatZeroRowsColumns_HYPRE(Mat A, PetscInt numRows, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
1456 {
1457   hypre_ParCSRMatrix *parcsr;
1458   HYPRE_Int          *lrows;
1459   PetscInt           rst,ren,i;
1460   PetscErrorCode     ierr;
1461 
1462   PetscFunctionBegin;
1463   if (x || b) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"To be implemented");
1464   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1465   ierr = PetscMalloc1(numRows,&lrows);CHKERRQ(ierr);
1466   ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr);
1467   for (i=0;i<numRows;i++) {
1468     if (rows[i] < rst || rows[i] >= ren)
1469       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Non-local rows not yet supported");
1470     lrows[i] = rows[i] - rst;
1471   }
1472   PetscStackCallStandard(hypre_ParCSRMatrixEliminateRowsCols,(parcsr,numRows,lrows));
1473   ierr = PetscFree(lrows);CHKERRQ(ierr);
1474   PetscFunctionReturn(0);
1475 }
1476 
1477 static PetscErrorCode MatZeroEntries_HYPRE_CSRMatrix(hypre_CSRMatrix *ha)
1478 {
1479   PetscErrorCode      ierr;
1480 
1481   PetscFunctionBegin;
1482   if (ha) {
1483     HYPRE_Int     *ii, size;
1484     HYPRE_Complex *a;
1485 
1486     size = hypre_CSRMatrixNumRows(ha);
1487     a    = hypre_CSRMatrixData(ha);
1488     ii   = hypre_CSRMatrixI(ha);
1489 
1490     if (a) { ierr = PetscMemzero(a,(ii[size])*sizeof(HYPRE_Complex));CHKERRQ(ierr); }
1491   }
1492   PetscFunctionReturn(0);
1493 }
1494 
1495 PetscErrorCode MatZeroEntries_HYPRE(Mat A)
1496 {
1497   hypre_ParCSRMatrix  *parcsr;
1498   PetscErrorCode      ierr;
1499 
1500   PetscFunctionBegin;
1501   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1502   /* diagonal part */
1503   ierr = MatZeroEntries_HYPRE_CSRMatrix(hypre_ParCSRMatrixDiag(parcsr));CHKERRQ(ierr);
1504   /* off-diagonal part */
1505   ierr = MatZeroEntries_HYPRE_CSRMatrix(hypre_ParCSRMatrixOffd(parcsr));CHKERRQ(ierr);
1506   PetscFunctionReturn(0);
1507 }
1508 
1509 
1510 static PetscErrorCode MatZeroRows_HYPRE_CSRMatrix(hypre_CSRMatrix *hA,PetscInt N,const PetscInt rows[],PetscScalar diag)
1511 {
1512   PetscInt        ii, jj, ibeg, iend, irow;
1513   PetscInt        *i, *j;
1514   PetscScalar     *a;
1515 
1516   PetscFunctionBegin;
1517 
1518   if (!hA) PetscFunctionReturn(0);
1519 
1520   i = hypre_CSRMatrixI(hA);
1521   j = hypre_CSRMatrixJ(hA);
1522   a = hypre_CSRMatrixData(hA);
1523 
1524   for (ii = 0; ii < N; ii++) {
1525     irow = rows[ii];
1526     ibeg = i[irow];
1527     iend = i[irow+1];
1528     for (jj = ibeg; jj < iend; jj++)
1529       if (j[jj] == irow) a[jj] = diag;
1530       else a[jj] = 0.0;
1531    }
1532 
1533    PetscFunctionReturn(0);
1534 }
1535 
1536 PetscErrorCode MatZeroRows_HYPRE(Mat A,PetscInt N,const PetscInt rows[],PetscScalar diag,Vec x,Vec b)
1537 {
1538   hypre_ParCSRMatrix  *parcsr;
1539   PetscInt            *lrows,len;
1540   PetscErrorCode      ierr;
1541 
1542   PetscFunctionBegin;
1543   if (x || b) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Does not support to modify the solution and the right hand size \n");
1544   /* retrieve the internal matrix */
1545   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1546   /* get locally owned rows */
1547   ierr = MatZeroRowsMapLocal_Private(A,N,rows,&len,&lrows);CHKERRQ(ierr);
1548   /* zero diagonal part */
1549   ierr = MatZeroRows_HYPRE_CSRMatrix(hypre_ParCSRMatrixDiag(parcsr),len,lrows,diag);CHKERRQ(ierr);
1550   /* zero off-diagonal part */
1551   ierr = MatZeroRows_HYPRE_CSRMatrix(hypre_ParCSRMatrixOffd(parcsr),len,lrows,0.0);CHKERRQ(ierr);
1552 
1553   ierr = PetscFree(lrows);CHKERRQ(ierr);
1554   PetscFunctionReturn(0);
1555 }
1556 
1557 
1558 PetscErrorCode MatAssemblyBegin_HYPRE(Mat mat,MatAssemblyType mode)
1559 {
1560   PetscErrorCode ierr;
1561 
1562   PetscFunctionBegin;
1563   if (mat->nooffprocentries) PetscFunctionReturn(0);
1564 
1565   ierr = MatStashScatterBegin_Private(mat,&mat->stash,mat->rmap->range);CHKERRQ(ierr);
1566   PetscFunctionReturn(0);
1567 }
1568 
1569 PetscErrorCode MatGetRow_HYPRE(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1570 {
1571   hypre_ParCSRMatrix  *parcsr;
1572   PetscErrorCode      ierr;
1573 
1574   PetscFunctionBegin;
1575   /* retrieve the internal matrix */
1576   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1577   /* call HYPRE API */
1578   PetscStackCallStandard(HYPRE_ParCSRMatrixGetRow,(parcsr,row,nz,idx,v));
1579   PetscFunctionReturn(0);
1580 }
1581 
1582 PetscErrorCode MatRestoreRow_HYPRE(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v)
1583 {
1584   hypre_ParCSRMatrix  *parcsr;
1585   PetscErrorCode      ierr;
1586 
1587   PetscFunctionBegin;
1588   /* retrieve the internal matrix */
1589   ierr = MatHYPREGetParCSR_HYPRE(A,&parcsr);CHKERRQ(ierr);
1590   /* call HYPRE API */
1591   PetscStackCallStandard(HYPRE_ParCSRMatrixRestoreRow,(parcsr,row,nz,idx,v));
1592   PetscFunctionReturn(0);
1593 }
1594 
1595 
1596 
1597 PetscErrorCode MatGetValues_HYPRE(Mat A,PetscInt m,const PetscInt idxm[],PetscInt n,const PetscInt idxn[],PetscScalar v[])
1598 {
1599   HYPRE_IJMatrix     *hIJ = (HYPRE_IJMatrix*)A->data;
1600   PetscErrorCode      ierr;
1601   PetscInt            i;
1602   PetscFunctionBegin;
1603   if (!m || !n) PetscFunctionReturn(0);
1604 
1605   /* Ignore negative row indices
1606    * And negative column indices should be automatically ignored in hypre
1607    * */
1608   for (i=0; i<m; i++)
1609     if (idxm[i] >= 0) PetscStackCallStandard(HYPRE_IJMatrixGetValues,(*hIJ,1,&n,(HYPRE_Int*)&idxm[i],(HYPRE_Int*)idxn,&v[i*n]));
1610 
1611   PetscFunctionReturn(0);
1612 }
1613 
1614 
1615 /*MC
1616    MATHYPRE - MATHYPRE = "hypre" - A matrix type to be used for sequential and parallel sparse matrices
1617           based on the hypre IJ interface.
1618 
1619    Level: intermediate
1620 
1621 .seealso: MatCreate()
1622 M*/
1623 
1624 PETSC_EXTERN PetscErrorCode MatCreate_HYPRE(Mat B)
1625 {
1626   Mat_HYPRE      *hB;
1627   PetscErrorCode ierr;
1628 
1629   PetscFunctionBegin;
1630   ierr           = PetscNewLog(B,&hB);CHKERRQ(ierr);
1631   hB->inner_free = PETSC_TRUE;
1632   hB->available  = PETSC_TRUE;
1633   hB->size       = 0;
1634   hB->array      = NULL;
1635 
1636   B->data       = (void*)hB;
1637   B->rmap->bs   = 1;
1638   B->assembled  = PETSC_FALSE;
1639 
1640   ierr = PetscMemzero(B->ops,sizeof(struct _MatOps));CHKERRQ(ierr);
1641   B->ops->mult            = MatMult_HYPRE;
1642   B->ops->multtranspose   = MatMultTranspose_HYPRE;
1643   B->ops->setup           = MatSetUp_HYPRE;
1644   B->ops->destroy         = MatDestroy_HYPRE;
1645 
1646   /* build cache for off array entries formed */
1647   ierr = MatStashCreate_Private(PetscObjectComm((PetscObject)B),1,&B->stash);CHKERRQ(ierr);
1648   B->ops->assemblyend     = MatAssemblyEnd_HYPRE;
1649   B->ops->assemblybegin   = MatAssemblyBegin_HYPRE;
1650 
1651   B->ops->ptap            = MatPtAP_HYPRE_HYPRE;
1652   B->ops->matmult         = MatMatMult_HYPRE_HYPRE;
1653   B->ops->setvalues       = MatSetValues_HYPRE;
1654   B->ops->missingdiagonal = MatMissingDiagonal_HYPRE;
1655   B->ops->scale           = MatScale_HYPRE;
1656   B->ops->zerorowscolumns = MatZeroRowsColumns_HYPRE;
1657   B->ops->zeroentries     = MatZeroEntries_HYPRE;
1658   B->ops->zerorows        = MatZeroRows_HYPRE;
1659   B->ops->getrow          = MatGetRow_HYPRE;
1660   B->ops->restorerow      = MatRestoreRow_HYPRE;
1661   B->ops->getvalues       = MatGetValues_HYPRE;
1662 
1663   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&hB->comm);CHKERRQ(ierr);
1664   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRE);CHKERRQ(ierr);
1665   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_aij_C",MatConvert_HYPRE_AIJ);CHKERRQ(ierr);
1666   ierr = PetscObjectComposeFunction((PetscObject)B,"MatConvert_hypre_is_C",MatConvert_HYPRE_IS);CHKERRQ(ierr);
1667   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_seqaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1668   ierr = PetscObjectComposeFunction((PetscObject)B,"MatPtAP_mpiaij_hypre_C",MatPtAP_AIJ_HYPRE);CHKERRQ(ierr);
1669   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPRESetPreallocation_C",MatHYPRESetPreallocation_HYPRE);CHKERRQ(ierr);
1670   ierr = PetscObjectComposeFunction((PetscObject)B,"MatHYPREGetParCSR_C",MatHYPREGetParCSR_HYPRE);CHKERRQ(ierr);
1671   PetscFunctionReturn(0);
1672 }
1673 
1674 static PetscErrorCode hypre_array_destroy(void *ptr)
1675 {
1676    PetscFunctionBegin;
1677    hypre_TFree(ptr,HYPRE_MEMORY_HOST);
1678    PetscFunctionReturn(0);
1679 }
1680