xref: /petsc/src/dm/impls/da/hypre/mhyp.c (revision 84df9cb40eca90ea9b18a456fab7a4ecc7f6c1a4) !
1 
2 /*
3     Creates hypre ijmatrix from PETSc matrix
4 */
5 #include <petscsys.h>
6 #include <private/matimpl.h>          /*I "petscmat.h" I*/
7 #include <petscdmda.h>                /*I "petscdmda.h" I*/
8 #include <../src/dm/impls/da/hypre/mhyp.h>
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "MatHYPRE_IJMatrixPreallocate"
12 PetscErrorCode MatHYPRE_IJMatrixPreallocate(Mat A_d, Mat A_o,HYPRE_IJMatrix ij)
13 {
14   PetscErrorCode ierr;
15   PetscInt       i;
16   PetscInt       n_d,*ia_d,n_o,*ia_o;
17   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
18   PetscInt       *nnz_d=PETSC_NULL,*nnz_o=PETSC_NULL;
19 
20   PetscFunctionBegin;
21   if (A_d) { /* determine number of nonzero entries in local diagonal part */
22     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,PETSC_NULL,&done_d);CHKERRQ(ierr);
23     if (done_d) {
24       ierr = PetscMalloc(n_d*sizeof(PetscInt),&nnz_d);CHKERRQ(ierr);
25       for (i=0; i<n_d; i++) {
26         nnz_d[i] = ia_d[i+1] - ia_d[i];
27       }
28     }
29     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,PETSC_NULL,&done_d);CHKERRQ(ierr);
30   }
31   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
32     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,PETSC_NULL,&done_o);CHKERRQ(ierr);
33     if (done_o) {
34       ierr = PetscMalloc(n_o*sizeof(PetscInt),&nnz_o);CHKERRQ(ierr);
35       for (i=0; i<n_o; i++) {
36         nnz_o[i] = ia_o[i+1] - ia_o[i];
37       }
38     }
39     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,PETSC_NULL,&done_o);CHKERRQ(ierr);
40   }
41   if (done_d) {    /* set number of nonzeros in HYPRE IJ matrix */
42     if (!done_o) { /* only diagonal part */
43       ierr = PetscMalloc(n_d*sizeof(PetscInt),&nnz_o);CHKERRQ(ierr);
44       for (i=0; i<n_d; i++) {
45         nnz_o[i] = 0;
46       }
47     }
48     PetscStackCallHypre(0,HYPRE_IJMatrixSetDiagOffdSizes,(ij,nnz_d,nnz_o));
49     ierr = PetscFree(nnz_d);CHKERRQ(ierr);
50     ierr = PetscFree(nnz_o);CHKERRQ(ierr);
51   }
52   PetscFunctionReturn(0);
53 }
54 
55 #undef __FUNCT__
56 #define __FUNCT__ "MatHYPRE_IJMatrixCreate"
57 PetscErrorCode MatHYPRE_IJMatrixCreate(Mat A,HYPRE_IJMatrix *ij)
58 {
59   PetscErrorCode ierr;
60   int            rstart,rend,cstart,cend;
61 
62   PetscFunctionBegin;
63   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
64   PetscValidType(A,1);
65   PetscValidPointer(ij,2);
66   ierr = MatPreallocated(A);CHKERRQ(ierr);
67   rstart = A->rmap->rstart;
68   rend   = A->rmap->rend;
69   cstart = A->cmap->rstart;
70   cend   = A->cmap->rend;
71   PetscStackCallHypre(0,HYPRE_IJMatrixCreate,(((PetscObject)A)->comm,rstart,rend-1,cstart,cend-1,ij));
72   PetscStackCallHypre(0,HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
73   {
74     PetscBool   same;
75     Mat         A_d,A_o;
76     PetscInt    *colmap;
77     ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&same);CHKERRQ(ierr);
78     if (same) {
79       ierr = MatMPIAIJGetSeqAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
80       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,*ij);CHKERRQ(ierr);
81       PetscFunctionReturn(0);
82     }
83     ierr = PetscTypeCompare((PetscObject)A,MATMPIBAIJ,&same);CHKERRQ(ierr);
84     if (same) {
85       ierr = MatMPIBAIJGetSeqBAIJ(A,&A_d,&A_o,&colmap);CHKERRQ(ierr);
86       ierr = MatHYPRE_IJMatrixPreallocate(A_d,A_o,*ij);CHKERRQ(ierr);
87       PetscFunctionReturn(0);
88     }
89     ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
90     if (same) {
91       ierr = MatHYPRE_IJMatrixPreallocate(A,PETSC_NULL,*ij);CHKERRQ(ierr);
92       PetscFunctionReturn(0);
93     }
94     ierr = PetscTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
95     if (same) {
96       ierr = MatHYPRE_IJMatrixPreallocate(A,PETSC_NULL,*ij);CHKERRQ(ierr);
97       PetscFunctionReturn(0);
98     }
99   }
100   PetscFunctionReturn(0);
101 }
102 
103 extern PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat,HYPRE_IJMatrix);
104 extern PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat,HYPRE_IJMatrix);
105 /*
106     Copies the data over (column indices, numerical values) to hypre matrix
107 */
108 
109 #undef __FUNCT__
110 #define __FUNCT__ "MatHYPRE_IJMatrixCopy"
111 PetscErrorCode MatHYPRE_IJMatrixCopy(Mat A,HYPRE_IJMatrix ij)
112 {
113   PetscErrorCode    ierr;
114   PetscInt          i,rstart,rend,ncols;
115   const PetscScalar *values;
116   const PetscInt    *cols;
117   PetscBool         flg;
118 
119   PetscFunctionBegin;
120   ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
121   if (flg) {
122     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
123     PetscFunctionReturn(0);
124   }
125   ierr = PetscTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
126   if (flg) {
127     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
128     PetscFunctionReturn(0);
129   }
130 
131   ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
132   PetscStackCallHypre(0,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     PetscStackCallHypre(0,HYPRE_IJMatrixSetValues,(ij,1,&ncols,&i,cols,values));
137     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
138   }
139   PetscStackCallHypre(0,HYPRE_IJMatrixAssemble,(ij));
140   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
141   PetscFunctionReturn(0);
142 }
143 
144 /*
145     This copies the CSR format directly from the PETSc data structure to the hypre
146     data structure without calls to MatGetRow() or hypre's set values.
147 
148 */
149 #include <_hypre_IJ_mv.h>
150 #include <HYPRE_IJ_mv.h>
151 #include <../src/mat/impls/aij/mpi/mpiaij.h>
152 
153 #undef __FUNCT__
154 #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ"
155 PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A,HYPRE_IJMatrix ij)
156 {
157   PetscErrorCode        ierr;
158   Mat_SeqAIJ            *pdiag = (Mat_SeqAIJ*)A->data;
159 
160   hypre_ParCSRMatrix    *par_matrix;
161   hypre_AuxParCSRMatrix *aux_matrix;
162   hypre_CSRMatrix       *hdiag;
163 
164   PetscFunctionBegin;
165   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
166   PetscValidType(A,1);
167   PetscValidPointer(ij,2);
168 
169   ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
170   PetscStackCallHypre(0,HYPRE_IJMatrixInitialize,(ij));
171   par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij);
172   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
173   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
174 
175   /*
176        this is the Hack part where we monkey directly with the hypre datastructures
177   */
178   ierr = PetscMemcpy(hdiag->i,pdiag->i,(A->rmap->n + 1)*sizeof(PetscInt));
179   ierr = PetscMemcpy(hdiag->j,pdiag->j,pdiag->nz*sizeof(PetscInt));
180   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
181 
182   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
183   PetscStackCallHypre(0,HYPRE_IJMatrixAssemble,(ij));
184   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
185   PetscFunctionReturn(0);
186 }
187 
188 #undef __FUNCT__
189 #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_MPIAIJ"
190 PetscErrorCode MatHYPRE_IJMatrixFastCopy_MPIAIJ(Mat A,HYPRE_IJMatrix ij)
191 {
192   PetscErrorCode        ierr;
193   Mat_MPIAIJ            *pA = (Mat_MPIAIJ*)A->data;
194   Mat_SeqAIJ            *pdiag,*poffd;
195   PetscInt              i,*garray = pA->garray,*jj,cstart,*pjj;
196 
197   hypre_ParCSRMatrix    *par_matrix;
198   hypre_AuxParCSRMatrix *aux_matrix;
199   hypre_CSRMatrix       *hdiag,*hoffd;
200 
201   PetscFunctionBegin;
202   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
203   PetscValidType(A,1);
204   PetscValidPointer(ij,2);
205   pdiag = (Mat_SeqAIJ*) pA->A->data;
206   poffd = (Mat_SeqAIJ*) pA->B->data;
207   /* cstart is only valid for square MPIAIJ layed out in the usual way */
208   ierr = MatGetOwnershipRange(A,&cstart,PETSC_NULL);CHKERRQ(ierr);
209 
210   ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
211   PetscStackCallHypre(0,HYPRE_IJMatrixInitialize,(ij));
212   par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij);
213   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
214   hdiag = hypre_ParCSRMatrixDiag(par_matrix);
215   hoffd = hypre_ParCSRMatrixOffd(par_matrix);
216 
217   /*
218        this is the Hack part where we monkey directly with the hypre datastructures
219   */
220   ierr = PetscMemcpy(hdiag->i,pdiag->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
221   /* need to shift the diag column indices (hdiag->j) back to global numbering since hypre is expecting this */
222   jj  = hdiag->j;
223   pjj = pdiag->j;
224   for (i=0; i<pdiag->nz; i++) {
225     jj[i] = cstart + pjj[i];
226   }
227   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
228   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
229   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
230      If we hacked a hypre a bit more we might be able to avoid this step */
231   jj  = hoffd->j;
232   pjj = poffd->j;
233   for (i=0; i<poffd->nz; i++) {
234     jj[i] = garray[pjj[i]];
235   }
236   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
237 
238   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
239   PetscStackCallHypre(0,HYPRE_IJMatrixAssemble,(ij));
240   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243 
244 /*
245     Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format
246 
247     This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first
248     which will corrupt the PETSc data structure if we did this. Need a work around to this problem.
249 */
250 #include <_hypre_IJ_mv.h>
251 #include <HYPRE_IJ_mv.h>
252 
253 #undef __FUNCT__
254 #define __FUNCT__ "MatHYPRE_IJMatrixLink"
255 PetscErrorCode MatHYPRE_IJMatrixLink(Mat A,HYPRE_IJMatrix *ij)
256 {
257   PetscErrorCode        ierr;
258   int                   rstart,rend,cstart,cend;
259   PetscBool             flg;
260   hypre_AuxParCSRMatrix *aux_matrix;
261 
262   PetscFunctionBegin;
263   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
264   PetscValidType(A,1);
265   PetscValidPointer(ij,2);
266   ierr = PetscTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
267   if (!flg) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices");
268   ierr = MatPreallocated(A);CHKERRQ(ierr);
269 
270   ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
271   rstart = A->rmap->rstart;
272   rend   = A->rmap->rend;
273   cstart = A->cmap->rstart;
274   cend   = A->cmap->rend;
275   PetscStackCallHypre(0,HYPRE_IJMatrixCreate,(((PetscObject)A)->comm,rstart,rend-1,cstart,cend-1,ij));
276   PetscStackCallHypre(0,HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
277 
278   PetscStackCallHypre(0,HYPRE_IJMatrixInitialize,(*ij));
279   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij);
280 
281   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
282 
283   /* this is the Hack part where we monkey directly with the hypre datastructures */
284 
285   PetscStackCallHypre(0,HYPRE_IJMatrixAssemble,(*ij));
286   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
287   PetscFunctionReturn(0);
288 }
289 
290 /* -----------------------------------------------------------------------------------------------------------------*/
291 
292 /*MC
293    MATHYPRESTRUCT - MATHYPRESTRUCT = "hyprestruct" - A matrix type to be used for parallel sparse matrices
294           based on the hypre HYPRE_StructMatrix.
295 
296    Level: intermediate
297 
298    Notes: Unlike the more general support for blocks in hypre this allows only one block per process and requires the block
299           be defined by a DMDA.
300 
301           The matrix needs a DMDA associated with it by either a call to MatSetDM() or if the matrix is obtained from DMGetMatrix()
302 
303 .seealso: MatCreate(), PCPFMG, MatSetDM(), DMGetMatrix()
304 M*/
305 
306 
307 #undef __FUNCT__
308 #define __FUNCT__ "MatSetValuesLocal_HYPREStruct_3d"
309 PetscErrorCode  MatSetValuesLocal_HYPREStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
310 {
311   PetscErrorCode    ierr;
312   PetscInt          i,j,stencil,index[3],row,entries[7];
313   const PetscScalar *values = y;
314   Mat_HYPREStruct   *ex = (Mat_HYPREStruct*) mat->data;
315 
316   PetscFunctionBegin;
317   for (i=0; i<nrow; i++) {
318     for (j=0; j<ncol; j++) {
319       stencil = icol[j] - irow[i];
320       if (!stencil) {
321         entries[j] = 3;
322       } else if (stencil == -1) {
323         entries[j] = 2;
324       } else if (stencil == 1) {
325         entries[j] = 4;
326       } else if (stencil == -ex->gnx) {
327         entries[j] = 1;
328       } else if (stencil == ex->gnx) {
329         entries[j] = 5;
330       } else if (stencil == -ex->gnxgny) {
331         entries[j] = 0;
332       } else if (stencil == ex->gnxgny) {
333         entries[j] = 6;
334       } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row %D local column %D have bad stencil %D",irow[i],icol[j],stencil);
335     }
336     row = ex->gindices[irow[i]] - ex->rstart;
337     index[0] = ex->xs + (row % ex->nx);
338     index[1] = ex->ys + ((row/ex->nx) % ex->ny);
339     index[2] = ex->zs + (row/(ex->nxny));
340     if (addv == ADD_VALUES) {
341       PetscStackCallHypre(0,HYPRE_StructMatrixAddToValues,(ex->hmat,index,ncol,entries,(PetscScalar*)values));
342     } else {
343       PetscStackCallHypre(0,HYPRE_StructMatrixSetValues,(ex->hmat,index,ncol,entries,(PetscScalar*)values));
344     }
345     values += ncol;
346   }
347   PetscFunctionReturn(0);
348 }
349 
350 #undef __FUNCT__
351 #define __FUNCT__ "MatZeroRowsLocal_HYPREStruct_3d"
352 PetscErrorCode  MatZeroRowsLocal_HYPREStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscScalar d,Vec x,Vec b)
353 {
354   PetscErrorCode    ierr;
355   PetscInt          i,index[3],row,entries[7] = {0,1,2,3,4,5,6};
356   PetscScalar       values[7];
357   Mat_HYPREStruct   *ex = (Mat_HYPREStruct*) mat->data;
358 
359   PetscFunctionBegin;
360   if (x && b) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"No support");
361   ierr = PetscMemzero(values,7*sizeof(PetscScalar));CHKERRQ(ierr);
362   values[3] = d;
363   for (i=0; i<nrow; i++) {
364     row = ex->gindices[irow[i]] - ex->rstart;
365     index[0] = ex->xs + (row % ex->nx);
366     index[1] = ex->ys + ((row/ex->nx) % ex->ny);
367     index[2] = ex->zs + (row/(ex->nxny));
368     PetscStackCallHypre(0,HYPRE_StructMatrixSetValues,(ex->hmat,index,7,entries,values));
369   }
370   PetscStackCallHypre(0,HYPRE_StructMatrixAssemble,(ex->hmat));
371   PetscFunctionReturn(0);
372 }
373 
374 #undef __FUNCT__
375 #define __FUNCT__ "MatZeroEntries_HYPREStruct_3d"
376 PetscErrorCode MatZeroEntries_HYPREStruct_3d(Mat mat)
377 {
378   PetscErrorCode ierr;
379   PetscInt       indices[7] = {0,1,2,3,4,5,6};
380   Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data;
381 
382   PetscFunctionBegin;
383   /* hypre has no public interface to do this */
384   PetscStackCallHypre(0,hypre_StructMatrixClearBoxValues,(ex->hmat,&ex->hbox,7,indices,0,1));
385   PetscStackCallHypre(0,HYPRE_StructMatrixAssemble,(ex->hmat));
386   PetscFunctionReturn(0);
387 }
388 
389 #undef __FUNCT__
390 #define __FUNCT__ "MatSetDM_HYPREStruct"
391 PetscErrorCode  MatSetDM_HYPREStruct(Mat mat,DM da)
392 {
393   PetscErrorCode  ierr;
394   Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data;
395   PetscInt         dim,dof,sw[3],nx,ny,nz;
396   int              ilower[3],iupper[3],ssize,i;
397   DMDABoundaryType   px,py,pz;
398   DMDAStencilType    st;
399 
400   PetscFunctionBegin;
401   ex->da = da;
402   ierr   = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);
403 
404   ierr = DMDAGetInfo(ex->da,&dim,0,0,0,0,0,0,&dof,&sw[0],&px,&py,&pz,&st);CHKERRQ(ierr);
405   ierr = DMDAGetCorners(ex->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
406   iupper[0] += ilower[0] - 1;
407   iupper[1] += ilower[1] - 1;
408   iupper[2] += ilower[2] - 1;
409 
410   /* the hypre_Box is used to zero out the matrix entries in MatZeroValues() */
411   ex->hbox.imin[0] = ilower[0];
412   ex->hbox.imin[1] = ilower[1];
413   ex->hbox.imin[2] = ilower[2];
414   ex->hbox.imax[0] = iupper[0];
415   ex->hbox.imax[1] = iupper[1];
416   ex->hbox.imax[2] = iupper[2];
417 
418   /* create the hypre grid object and set its information */
419   if (dof > 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Currently only support for scalar problems");
420   if (px || py || pz) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add periodic support by calling HYPRE_StructGridSetPeriodic()");
421   PetscStackCallHypre(0,HYPRE_StructGridCreate,(ex->hcomm,dim,&ex->hgrid));
422   PetscStackCallHypre(0,HYPRE_StructGridSetExtents,(ex->hgrid,ilower,iupper));
423   PetscStackCallHypre(0,HYPRE_StructGridAssemble,(ex->hgrid));
424 
425   sw[1] = sw[0];
426   sw[2] = sw[1];
427   PetscStackCallHypre(0,HYPRE_StructGridSetNumGhost,(ex->hgrid,sw));
428 
429   /* create the hypre stencil object and set its information */
430   if (sw[0] > 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for wider stencils");
431   if (st == DMDA_STENCIL_BOX) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for box stencils");
432   if (dim == 1) {
433     int offsets[3][1] = {{-1},{0},{1}};
434     ssize = 3;
435     PetscStackCallHypre(0,HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil));
436     for (i=0; i<ssize; i++) {
437       PetscStackCallHypre(0,HYPRE_StructStencilSetElement,(ex->hstencil,i,offsets[i]));
438     }
439   } else if (dim == 2) {
440     int offsets[5][2] = {{0,-1},{-1,0},{0,0},{1,0},{0,1}};
441     ssize = 5;
442     PetscStackCallHypre(0,HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil));
443     for (i=0; i<ssize; i++) {
444       PetscStackCallHypre(0,HYPRE_StructStencilSetElement,(ex->hstencil,i,offsets[i]));
445     }
446   } else if (dim == 3) {
447     int offsets[7][3] = {{0,0,-1},{0,-1,0},{-1,0,0},{0,0,0},{1,0,0},{0,1,0},{0,0,1}};
448     ssize = 7;
449     PetscStackCallHypre(0,HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil));
450     for (i=0; i<ssize; i++) {
451       PetscStackCallHypre(0,HYPRE_StructStencilSetElement,(ex->hstencil,i,offsets[i]));
452     }
453   }
454 
455   /* create the HYPRE vector for rhs and solution */
456   PetscStackCallHypre(0,HYPRE_StructVectorCreate,(ex->hcomm,ex->hgrid,&ex->hb));
457   PetscStackCallHypre(0,HYPRE_StructVectorCreate,(ex->hcomm,ex->hgrid,&ex->hx));
458   PetscStackCallHypre(0,HYPRE_StructVectorInitialize,(ex->hb));
459   PetscStackCallHypre(0,HYPRE_StructVectorInitialize,(ex->hx));
460   PetscStackCallHypre(0,HYPRE_StructVectorAssemble,(ex->hb));
461   PetscStackCallHypre(0,HYPRE_StructVectorAssemble,(ex->hx));
462 
463   /* create the hypre matrix object and set its information */
464   PetscStackCallHypre(0,HYPRE_StructMatrixCreate,(ex->hcomm,ex->hgrid,ex->hstencil,&ex->hmat));
465   PetscStackCallHypre(0,HYPRE_StructGridDestroy,(ex->hgrid));
466   PetscStackCallHypre(0,HYPRE_StructStencilDestroy,(ex->hstencil));
467   if (ex->needsinitialization) {
468     PetscStackCallHypre(0,HYPRE_StructMatrixInitialize,(ex->hmat));
469     ex->needsinitialization = PETSC_FALSE;
470   }
471 
472   /* set the global and local sizes of the matrix */
473   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
474   ierr = MatSetSizes(mat,dof*nx*ny*nz,dof*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
475   ierr = PetscLayoutSetBlockSize(mat->rmap,1);CHKERRQ(ierr);
476   ierr = PetscLayoutSetBlockSize(mat->cmap,1);CHKERRQ(ierr);
477   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
478   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
479 
480   if (dim == 3) {
481     mat->ops->setvalueslocal = MatSetValuesLocal_HYPREStruct_3d;
482     mat->ops->zerorowslocal  = MatZeroRowsLocal_HYPREStruct_3d;
483     mat->ops->zeroentries    = MatZeroEntries_HYPREStruct_3d;
484     ierr = MatZeroEntries_HYPREStruct_3d(mat);CHKERRQ(ierr);
485   } else SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only support for 3d DMDA currently");
486 
487   /* get values that will be used repeatedly in MatSetValuesLocal() and MatZeroRowsLocal() repeatedly */
488   ierr = MatGetOwnershipRange(mat,&ex->rstart,PETSC_NULL);CHKERRQ(ierr);
489   ierr = DMDAGetGlobalIndices(ex->da,PETSC_NULL,&ex->gindices);CHKERRQ(ierr);
490   ierr = DMDAGetGhostCorners(ex->da,0,0,0,&ex->gnx,&ex->gnxgny,0);CHKERRQ(ierr);
491   ex->gnxgny *= ex->gnx;
492   ierr = DMDAGetCorners(ex->da,&ex->xs,&ex->ys,&ex->zs,&ex->nx,&ex->ny,0);CHKERRQ(ierr);
493   ex->nxny = ex->nx*ex->ny;
494   PetscFunctionReturn(0);
495 }
496 
497 #undef __FUNCT__
498 #define __FUNCT__ "MatMult_HYPREStruct"
499 PetscErrorCode MatMult_HYPREStruct(Mat A,Vec x,Vec y)
500 {
501   PetscErrorCode  ierr;
502   PetscScalar     *xx,*yy;
503   int             ilower[3],iupper[3];
504   Mat_HYPREStruct *mx = (Mat_HYPREStruct *)(A->data);
505 
506   PetscFunctionBegin;
507   ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
508   iupper[0] += ilower[0] - 1;
509   iupper[1] += ilower[1] - 1;
510   iupper[2] += ilower[2] - 1;
511 
512   /* copy x values over to hypre */
513   PetscStackCallHypre(0,HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
514   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
515   PetscStackCallHypre(0,HYPRE_StructVectorSetBoxValues,(mx->hb,ilower,iupper,xx));
516   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
517   PetscStackCallHypre(0,HYPRE_StructVectorAssemble,(mx->hb));
518   PetscStackCallHypre(0,HYPRE_StructMatrixMatvec,(1.0,mx->hmat,mx->hb,0.0,mx->hx));
519 
520   /* copy solution values back to PETSc */
521   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
522   PetscStackCallHypre(0,HYPRE_StructVectorGetBoxValues,(mx->hx,ilower,iupper,yy));
523   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
524   PetscFunctionReturn(0);
525 }
526 
527 #undef __FUNCT__
528 #define __FUNCT__ "MatAssemblyEnd_HYPREStruct"
529 PetscErrorCode MatAssemblyEnd_HYPREStruct(Mat mat,MatAssemblyType mode)
530 {
531   Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data;
532   PetscErrorCode   ierr;
533 
534   PetscFunctionBegin;
535   PetscStackCallHypre(0,HYPRE_StructMatrixAssemble,(ex->hmat));
536   /* PetscStackCallHypre(0,HYPRE_StructMatrixPrint,("dummy",ex->hmat,0)); */
537   PetscFunctionReturn(0);
538 }
539 
540 #undef __FUNCT__
541 #define __FUNCT__ "MatZeroEntries_HYPREStruct"
542 PetscErrorCode MatZeroEntries_HYPREStruct(Mat mat)
543 {
544   PetscFunctionBegin;
545   /* before the DMDA is set to the matrix the zero doesn't need to do anything */
546   PetscFunctionReturn(0);
547 }
548 
549 
550 #undef __FUNCT__
551 #define __FUNCT__ "MatDestroy_HYPREStruct"
552 PetscErrorCode MatDestroy_HYPREStruct(Mat mat)
553 {
554   Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data;
555   PetscErrorCode  ierr;
556 
557   PetscFunctionBegin;
558   PetscStackCallHypre(0,HYPRE_StructMatrixDestroy,(ex->hmat));
559   PetscStackCallHypre(0,HYPRE_StructVectorDestroy,(ex->hx));
560   PetscStackCallHypre(0,HYPRE_StructVectorDestroy,(ex->hb));
561   PetscFunctionReturn(0);
562 }
563 
564 
565 EXTERN_C_BEGIN
566 #undef __FUNCT__
567 #define __FUNCT__ "MatCreate_HYPREStruct"
568 PetscErrorCode  MatCreate_HYPREStruct(Mat B)
569 {
570   Mat_HYPREStruct *ex;
571   PetscErrorCode  ierr;
572 
573   PetscFunctionBegin;
574   ierr            = PetscNewLog(B,Mat_HYPREStruct,&ex);CHKERRQ(ierr);
575   B->data         = (void*)ex;
576   B->rmap->bs     = 1;
577   B->assembled    = PETSC_FALSE;
578 
579   B->insertmode   = NOT_SET_VALUES;
580 
581   B->ops->assemblyend    = MatAssemblyEnd_HYPREStruct;
582   B->ops->mult           = MatMult_HYPREStruct;
583   B->ops->zeroentries    = MatZeroEntries_HYPREStruct;
584   B->ops->destroy        = MatDestroy_HYPREStruct;
585 
586   ex->needsinitialization = PETSC_TRUE;
587 
588   ierr = MPI_Comm_dup(((PetscObject)B)->comm,&(ex->hcomm));CHKERRQ(ierr);
589   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSetDM_C","MatSetDM_HYPREStruct",MatSetDM_HYPREStruct);CHKERRQ(ierr);
590   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRESTRUCT);CHKERRQ(ierr);
591   PetscFunctionReturn(0);
592 }
593 EXTERN_C_END
594 
595 
596 /*MC
597    MATHYPRESSTRUCT - MATHYPRESSTRUCT = "hypresstruct" - A matrix type to be used for parallel sparse matrices
598           based on the hypre HYPRE_SStructMatrix.
599 
600 
601    Level: intermediate
602 
603    Notes: Unlike hypre's general semi-struct object consisting of a collection of structured-grid objects and unstructured
604           grid objects, we will restrict the semi-struct objects to consist of only structured-grid components.
605 
606           Unlike the more general support for parts and blocks in hypre this allows only one part, and one block per process and requires the block
607           be defined by a DMDA.
608 
609           The matrix needs a DMDA associated with it by either a call to MatSetDM() or if the matrix is obtained from DMGetMatrix()
610 
611 M*/
612 
613 #undef __FUNCT__
614 #define __FUNCT__ "MatSetValuesLocal_HYPRESStruct_3d"
615 PetscErrorCode  MatSetValuesLocal_HYPRESStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
616 {
617   PetscErrorCode    ierr;
618   PetscInt          i,j,stencil,index[3];
619   const PetscScalar *values = y;
620   Mat_HYPRESStruct  *ex = (Mat_HYPRESStruct*) mat->data;
621 
622   int               part= 0; /* Petsc sstruct interface only allows 1 part */
623   int               ordering;
624   int               grid_rank, to_grid_rank;
625   int               var_type, to_var_type;
626   int               to_var_entry = 0;
627 
628   int               nvars= ex->nvars;
629   PetscInt          row,*entries;
630 
631   PetscFunctionBegin;
632   ierr = PetscMalloc(7*nvars*sizeof(PetscInt),&entries);CHKERRQ(ierr);
633   ordering= ex-> dofs_order; /* ordering= 0   nodal ordering
634                                           1   variable ordering */
635   /* stencil entries are orderer by variables: var0_stencil0, var0_stencil1, ..., var0_stencil6, var1_stencil0, var1_stencil1, ...  */
636 
637   if (!ordering) {
638     for (i=0; i<nrow; i++) {
639       grid_rank= irow[i]/nvars;
640       var_type = (irow[i] % nvars);
641 
642       for (j=0; j<ncol; j++) {
643         to_grid_rank= icol[j]/nvars;
644         to_var_type = (icol[j] % nvars);
645 
646         to_var_entry= to_var_entry*7;
647         entries[j]= to_var_entry;
648 
649         stencil = to_grid_rank-grid_rank;
650         if (!stencil) {
651           entries[j] += 3;
652         } else if (stencil == -1) {
653           entries[j] += 2;
654         } else if (stencil == 1) {
655           entries[j] += 4;
656         } else if (stencil == -ex->gnx) {
657           entries[j] += 1;
658         } else if (stencil == ex->gnx) {
659           entries[j] += 5;
660         } else if (stencil == -ex->gnxgny) {
661           entries[j] += 0;
662         } else if (stencil == ex->gnxgny) {
663           entries[j] += 6;
664         } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row %D local column %D have bad stencil %D",irow[i],icol[j],stencil);
665       }
666 
667       row = ex->gindices[grid_rank] - ex->rstart;
668       index[0] = ex->xs + (row % ex->nx);
669       index[1] = ex->ys + ((row/ex->nx) % ex->ny);
670       index[2] = ex->zs + (row/(ex->nxny));
671 
672       if (addv == ADD_VALUES) {
673         PetscStackCallHypre(0,HYPRE_SStructMatrixAddToValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values));
674       } else {
675         PetscStackCallHypre(0,HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values));
676       }
677       values += ncol;
678     }
679   } else {
680     for (i=0; i<nrow; i++) {
681       var_type = irow[i]/(ex->gnxgnygnz);
682       grid_rank= irow[i] - var_type*(ex->gnxgnygnz);
683 
684       for (j=0; j<ncol; j++) {
685         to_var_type = icol[j]/(ex->gnxgnygnz);
686         to_grid_rank= icol[j] - to_var_type*(ex->gnxgnygnz);
687 
688         to_var_entry= to_var_entry*7;
689         entries[j]= to_var_entry;
690 
691         stencil = to_grid_rank-grid_rank;
692         if (!stencil) {
693           entries[j] += 3;
694         } else if (stencil == -1) {
695           entries[j] += 2;
696         } else if (stencil == 1) {
697           entries[j] += 4;
698         } else if (stencil == -ex->gnx) {
699           entries[j] += 1;
700         } else if (stencil == ex->gnx) {
701           entries[j] += 5;
702         } else if (stencil == -ex->gnxgny) {
703           entries[j] += 0;
704         } else if (stencil == ex->gnxgny) {
705           entries[j] += 6;
706         } else SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Local row %D local column %D have bad stencil %D",irow[i],icol[j],stencil);
707       }
708 
709       row = ex->gindices[grid_rank] - ex->rstart;
710       index[0] = ex->xs + (row % ex->nx);
711       index[1] = ex->ys + ((row/ex->nx) % ex->ny);
712       index[2] = ex->zs + (row/(ex->nxny));
713 
714       if (addv == ADD_VALUES) {
715         PetscStackCallHypre(0,HYPRE_SStructMatrixAddToValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values));
716       } else {
717         PetscStackCallHypre(0,HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,ncol,entries,(PetscScalar*)values));
718       }
719       values += ncol;
720     }
721   }
722   ierr = PetscFree(entries);CHKERRQ(ierr);
723   PetscFunctionReturn(0);
724 }
725 
726 #undef __FUNCT__
727 #define __FUNCT__ "MatZeroRowsLocal_HYPRESStruct_3d"
728 PetscErrorCode  MatZeroRowsLocal_HYPRESStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscScalar d,Vec x,Vec b)
729 {
730   PetscErrorCode    ierr;
731   PetscInt          i,index[3];
732   PetscScalar     **values;
733   Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data;
734 
735   int               part= 0; /* Petsc sstruct interface only allows 1 part */
736   int               ordering= ex->dofs_order;
737   int               grid_rank;
738   int               var_type;
739   int               nvars= ex->nvars;
740   PetscInt          row,*entries;
741 
742   PetscFunctionBegin;
743   if (x && b) SETERRQ(((PetscObject)mat)->comm,PETSC_ERR_SUP,"No support");
744   ierr = PetscMalloc(7*nvars*sizeof(PetscInt),&entries);CHKERRQ(ierr);
745 
746   ierr = PetscMalloc(nvars*sizeof(PetscScalar *),&values);CHKERRQ(ierr);
747   ierr = PetscMalloc(7*nvars*nvars*sizeof(PetscScalar),&values[0]);CHKERRQ(ierr);
748   for (i=1; i<nvars; i++) {
749      values[i] = values[i-1] + nvars*7;
750   }
751 
752   for (i=0; i< nvars; i++) {
753     ierr = PetscMemzero(values[i],nvars*7*sizeof(PetscScalar));CHKERRQ(ierr);
754     *(values[i]+3)= d;
755   }
756 
757   for (i= 0; i< nvars*7; i++) {
758     entries[i]= i;
759   }
760 
761   if (!ordering) {
762     for (i=0; i<nrow; i++) {
763        grid_rank= irow[i]/nvars;
764        var_type = (irow[i] % nvars);
765 
766        row = ex->gindices[grid_rank] - ex->rstart;
767        index[0] = ex->xs + (row % ex->nx);
768        index[1] = ex->ys + ((row/ex->nx) % ex->ny);
769        index[2] = ex->zs + (row/(ex->nxny));
770        PetscStackCallHypre(0,HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,7*nvars,entries,values[var_type]));
771     }
772   } else {
773     for (i=0; i<nrow; i++) {
774        var_type = irow[i]/(ex->gnxgnygnz);
775        grid_rank= irow[i] - var_type*(ex->gnxgnygnz);
776 
777        row = ex->gindices[grid_rank] - ex->rstart;
778        index[0] = ex->xs + (row % ex->nx);
779        index[1] = ex->ys + ((row/ex->nx) % ex->ny);
780        index[2] = ex->zs + (row/(ex->nxny));
781        PetscStackCallHypre(0,HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,index,var_type,7*nvars,entries,values[var_type]));
782     }
783   }
784   PetscStackCallHypre(0,HYPRE_SStructMatrixAssemble,(ex->ss_mat));
785   ierr = PetscFree(values[0]);CHKERRQ(ierr);
786   ierr = PetscFree(values);CHKERRQ(ierr);
787   ierr = PetscFree(entries);CHKERRQ(ierr);
788   PetscFunctionReturn(0);
789 }
790 
791 #undef __FUNCT__
792 #define __FUNCT__ "MatZeroEntries_HYPRESStruct_3d"
793 PetscErrorCode MatZeroEntries_HYPRESStruct_3d(Mat mat)
794 {
795   PetscErrorCode     ierr;
796   Mat_HYPRESStruct  *ex = (Mat_HYPRESStruct*) mat->data;
797   int                nvars= ex->nvars;
798   int                size;
799   int                part= 0; /* only one part */
800 
801   PetscFunctionBegin;
802   size = ((ex->hbox.imax[0])-(ex->hbox.imin[0])+1)*((ex->hbox.imax[1])-(ex->hbox.imin[1])+1)*((ex->hbox.imax[2])-(ex->hbox.imin[2])+1);
803   {
804      PetscInt          i,*entries;
805      PetscScalar      *values;
806      int               iupper[3], ilower[3];
807 
808      for (i= 0; i< 3; i++) {
809         ilower[i]= ex->hbox.imin[i];
810         iupper[i]= ex->hbox.imax[i];
811      }
812 
813      ierr = PetscMalloc2(nvars*7,PetscInt,&entries,nvars*7*size,PetscScalar,&values);CHKERRQ(ierr);
814      for (i= 0; i< nvars*7; i++) {
815         entries[i]= i;
816      }
817      ierr = PetscMemzero(values,nvars*7*size*sizeof(PetscScalar));CHKERRQ(ierr);
818 
819      for (i= 0; i< nvars; i++) {
820        PetscStackCallHypre(0,HYPRE_SStructMatrixSetBoxValues,(ex->ss_mat,part,ilower,iupper,i,nvars*7,entries,values));
821      }
822      ierr = PetscFree2(entries,values);CHKERRQ(ierr);
823   }
824   PetscStackCallHypre(0,HYPRE_SStructMatrixAssemble,(ex->ss_mat));
825   PetscFunctionReturn(0);
826 }
827 
828 
829 #undef __FUNCT__
830 #define __FUNCT__ "MatSetDM_HYPRESStruct"
831 PetscErrorCode  MatSetDM_HYPRESStruct(Mat mat,DM da)
832 {
833   PetscErrorCode    ierr;
834   Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data;
835   PetscInt          dim,dof,sw[3],nx,ny,nz;
836   int               ilower[3],iupper[3],ssize,i;
837   DMDABoundaryType    px,py,pz;
838   DMDAStencilType     st;
839   int               nparts= 1; /* assuming only one part */
840   int               part  = 0;
841 
842   PetscFunctionBegin;
843   ex->da = da;
844   ierr   = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);
845 
846   ierr = DMDAGetInfo(ex->da,&dim,0,0,0,0,0,0,&dof,&sw[0],&px,&py,&pz,&st);CHKERRQ(ierr);
847   ierr = DMDAGetCorners(ex->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
848   iupper[0] += ilower[0] - 1;
849   iupper[1] += ilower[1] - 1;
850   iupper[2] += ilower[2] - 1;
851   /* the hypre_Box is used to zero out the matrix entries in MatZeroValues() */
852   ex->hbox.imin[0] = ilower[0];
853   ex->hbox.imin[1] = ilower[1];
854   ex->hbox.imin[2] = ilower[2];
855   ex->hbox.imax[0] = iupper[0];
856   ex->hbox.imax[1] = iupper[1];
857   ex->hbox.imax[2] = iupper[2];
858 
859   ex->dofs_order   = 0;
860 
861   /* assuming that the same number of dofs on each gridpoint. Also assume all cell-centred based */
862   ex->nvars= dof;
863 
864   /* create the hypre grid object and set its information */
865   if (px || py || pz) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add periodic support by calling HYPRE_SStructGridSetPeriodic()");
866   PetscStackCallHypre(0,HYPRE_SStructGridCreate,(ex->hcomm,dim,nparts,&ex->ss_grid));
867 
868   PetscStackCallHypre(0,HYPRE_SStructGridSetExtents,(ex->ss_grid,part,ex->hbox.imin,ex->hbox.imax));
869 
870   {
871     HYPRE_SStructVariable *vartypes;
872     ierr = PetscMalloc(ex->nvars*sizeof(HYPRE_SStructVariable),&vartypes);CHKERRQ(ierr);
873     for (i= 0; i< ex->nvars; i++) {
874       vartypes[i]= HYPRE_SSTRUCT_VARIABLE_CELL;
875     }
876     PetscStackCallHypre(0,HYPRE_SStructGridSetVariables,(ex->ss_grid, part, ex->nvars,vartypes));
877     ierr = PetscFree(vartypes);CHKERRQ(ierr);
878   }
879   PetscStackCallHypre(0,HYPRE_SStructGridAssemble,(ex->ss_grid));
880 
881   sw[1] = sw[0];
882   sw[2] = sw[1];
883   /* PetscStackCallHypre(0,HYPRE_SStructGridSetNumGhost,(ex->ss_grid,sw)); */
884 
885   /* create the hypre stencil object and set its information */
886   if (sw[0] > 1) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for wider stencils");
887   if (st == DMDA_STENCIL_BOX) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Ask us to add support for box stencils");
888 
889   if (dim == 1) {
890     int offsets[3][1] = {{-1},{0},{1}};
891     int j, cnt;
892 
893     ssize = 3*(ex->nvars);
894     PetscStackCallHypre(0,HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil));
895     cnt= 0;
896     for (i= 0; i< (ex->nvars); i++) {
897        for (j= 0; j< 3; j++) {
898          PetscStackCallHypre(0,HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, offsets[j], i));
899           cnt++;
900        }
901     }
902   } else if (dim == 2) {
903     int offsets[5][2] = {{0,-1},{-1,0},{0,0},{1,0},{0,1}};
904     int j, cnt;
905 
906     ssize = 5*(ex->nvars);
907     PetscStackCallHypre(0,HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil));
908     cnt= 0;
909     for (i= 0; i< (ex->nvars); i++) {
910        for (j= 0; j< 5; j++) {
911          PetscStackCallHypre(0,HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, offsets[j], i));
912           cnt++;
913        }
914     }
915   } else if (dim == 3) {
916     int offsets[7][3] = {{0,0,-1},{0,-1,0},{-1,0,0},{0,0,0},{1,0,0},{0,1,0},{0,0,1}};
917     int j, cnt;
918 
919     ssize = 7*(ex->nvars);
920     PetscStackCallHypre(0,HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil));
921     cnt= 0;
922     for (i= 0; i< (ex->nvars); i++) {
923        for (j= 0; j< 7; j++) {
924          PetscStackCallHypre(0,HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, offsets[j], i));
925           cnt++;
926        }
927     }
928   }
929 
930   /* create the HYPRE graph */
931   PetscStackCallHypre(0,HYPRE_SStructGraphCreate,(ex->hcomm, ex->ss_grid, &(ex->ss_graph)));
932 
933   /* set the stencil graph. Note that each variable has the same graph. This means that each
934      variable couples to all the other variable and with the same stencil pattern. */
935   for (i= 0; i< (ex->nvars); i++) {
936     PetscStackCallHypre(0,HYPRE_SStructGraphSetStencil,(ex->ss_graph,part,i,ex->ss_stencil));
937   }
938   PetscStackCallHypre(0,HYPRE_SStructGraphAssemble,(ex->ss_graph));
939 
940   /* create the HYPRE sstruct vectors for rhs and solution */
941   PetscStackCallHypre(0,HYPRE_SStructVectorCreate,(ex->hcomm,ex->ss_grid,&ex->ss_b));
942   PetscStackCallHypre(0,HYPRE_SStructVectorCreate,(ex->hcomm,ex->ss_grid,&ex->ss_x));
943   PetscStackCallHypre(0,HYPRE_SStructVectorInitialize,(ex->ss_b));
944   PetscStackCallHypre(0,HYPRE_SStructVectorInitialize,(ex->ss_x));
945   PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(ex->ss_b));
946   PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(ex->ss_x));
947 
948   /* create the hypre matrix object and set its information */
949   PetscStackCallHypre(0,HYPRE_SStructMatrixCreate,(ex->hcomm,ex->ss_graph,&ex->ss_mat));
950   PetscStackCallHypre(0,HYPRE_SStructGridDestroy,(ex->ss_grid));
951   PetscStackCallHypre(0,HYPRE_SStructStencilDestroy,(ex->ss_stencil));
952   if (ex->needsinitialization) {
953     PetscStackCallHypre(0,HYPRE_SStructMatrixInitialize,(ex->ss_mat));
954     ex->needsinitialization = PETSC_FALSE;
955   }
956 
957   /* set the global and local sizes of the matrix */
958   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
959   ierr = MatSetSizes(mat,dof*nx*ny*nz,dof*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
960   ierr = PetscLayoutSetBlockSize(mat->rmap,1);CHKERRQ(ierr);
961   ierr = PetscLayoutSetBlockSize(mat->cmap,1);CHKERRQ(ierr);
962   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
963   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
964 
965   if (dim == 3) {
966     mat->ops->setvalueslocal = MatSetValuesLocal_HYPRESStruct_3d;
967     mat->ops->zerorowslocal  = MatZeroRowsLocal_HYPRESStruct_3d;
968     mat->ops->zeroentries    = MatZeroEntries_HYPRESStruct_3d;
969     ierr = MatZeroEntries_HYPRESStruct_3d(mat);CHKERRQ(ierr);
970   } else SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Only support for 3d DMDA currently");
971 
972   /* get values that will be used repeatedly in MatSetValuesLocal() and MatZeroRowsLocal() repeatedly */
973   ierr = MatGetOwnershipRange(mat,&ex->rstart,PETSC_NULL);CHKERRQ(ierr);
974   ierr = DMDAGetGlobalIndices(ex->da,PETSC_NULL,&ex->gindices);CHKERRQ(ierr);
975   ierr = DMDAGetGhostCorners(ex->da,0,0,0,&ex->gnx,&ex->gnxgny,&ex->gnxgnygnz);CHKERRQ(ierr);
976   ex->gnxgny    *= ex->gnx;
977   ex->gnxgnygnz *= ex->gnxgny;
978   ierr = DMDAGetCorners(ex->da,&ex->xs,&ex->ys,&ex->zs,&ex->nx,&ex->ny,&ex->nz);CHKERRQ(ierr);
979   ex->nxny   = ex->nx*ex->ny;
980   ex->nxnynz = ex->nz*ex->nxny;
981   PetscFunctionReturn(0);
982 }
983 
984 #undef __FUNCT__
985 #define __FUNCT__ "MatMult_HYPRESStruct"
986 PetscErrorCode MatMult_HYPRESStruct(Mat A,Vec x,Vec y)
987 {
988   PetscErrorCode    ierr;
989   PetscScalar      *xx,*yy;
990   int               ilower[3],iupper[3];
991   Mat_HYPRESStruct *mx = (Mat_HYPRESStruct *)(A->data);
992   int               ordering= mx->dofs_order;
993   int               nvars= mx->nvars;
994   int               part= 0;
995   int               size;
996   int               i;
997 
998   PetscFunctionBegin;
999   ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
1000   iupper[0] += ilower[0] - 1;
1001   iupper[1] += ilower[1] - 1;
1002   iupper[2] += ilower[2] - 1;
1003 
1004   size= 1;
1005   for (i= 0; i< 3; i++) {
1006      size*= (iupper[i]-ilower[i]+1);
1007   }
1008 
1009   /* copy x values over to hypre for variable ordering */
1010   if (ordering) {
1011     PetscStackCallHypre(0,HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1012      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1013      for (i= 0; i< nvars; i++) {
1014        PetscStackCallHypre(0,HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,xx+(size*i)));
1015      }
1016      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1017      PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(mx->ss_b));
1018      PetscStackCallHypre(0,HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
1019 
1020      /* copy solution values back to PETSc */
1021      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1022      for (i= 0; i< nvars; i++) {
1023        PetscStackCallHypre(0,HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,yy+(size*i)));
1024      }
1025      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1026   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
1027      PetscScalar     *z;
1028      int              j, k;
1029 
1030      ierr = PetscMalloc(nvars*size*sizeof(PetscScalar),&z);CHKERRQ(ierr);
1031      PetscStackCallHypre(0,HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1032      ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1033 
1034      /* transform nodal to hypre's variable ordering for sys_pfmg */
1035      for (i= 0; i< size; i++) {
1036         k= i*nvars;
1037         for (j= 0; j< nvars; j++) {
1038            z[j*size+i]= xx[k+j];
1039         }
1040      }
1041      for (i= 0; i< nvars; i++) {
1042        PetscStackCallHypre(0,HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,ilower,iupper,i,z+(size*i)));
1043      }
1044      ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1045      PetscStackCallHypre(0,HYPRE_SStructVectorAssemble,(mx->ss_b));
1046      PetscStackCallHypre(0,HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
1047 
1048      /* copy solution values back to PETSc */
1049      ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1050      for (i= 0; i< nvars; i++) {
1051        PetscStackCallHypre(0,HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,ilower,iupper,i,z+(size*i)));
1052      }
1053      /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
1054      for (i= 0; i< size; i++) {
1055         k= i*nvars;
1056         for (j= 0; j< nvars; j++) {
1057            yy[k+j]= z[j*size+i];
1058         }
1059      }
1060      ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1061      ierr = PetscFree(z);CHKERRQ(ierr);
1062   }
1063   PetscFunctionReturn(0);
1064 }
1065 
1066 #undef __FUNCT__
1067 #define __FUNCT__ "MatAssemblyEnd_HYPRESStruct"
1068 PetscErrorCode MatAssemblyEnd_HYPRESStruct(Mat mat,MatAssemblyType mode)
1069 {
1070   Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data;
1071   PetscErrorCode   ierr;
1072 
1073   PetscFunctionBegin;
1074   PetscStackCallHypre(0,HYPRE_SStructMatrixAssemble,(ex->ss_mat));
1075   PetscFunctionReturn(0);
1076 }
1077 
1078 #undef __FUNCT__
1079 #define __FUNCT__ "MatZeroEntries_HYPRESStruct"
1080 PetscErrorCode MatZeroEntries_HYPRESStruct(Mat mat)
1081 {
1082   PetscFunctionBegin;
1083   /* before the DMDA is set to the matrix the zero doesn't need to do anything */
1084   PetscFunctionReturn(0);
1085 }
1086 
1087 
1088 #undef __FUNCT__
1089 #define __FUNCT__ "MatDestroy_HYPRESStruct"
1090 PetscErrorCode MatDestroy_HYPRESStruct(Mat mat)
1091 {
1092   Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data;
1093   PetscErrorCode  ierr;
1094 
1095   PetscFunctionBegin;
1096   PetscStackCallHypre(0,HYPRE_SStructGraphDestroy,(ex->ss_graph));
1097   PetscStackCallHypre(0,HYPRE_SStructMatrixDestroy,(ex->ss_mat));
1098   PetscStackCallHypre(0,HYPRE_SStructVectorDestroy,(ex->ss_x));
1099   PetscStackCallHypre(0,HYPRE_SStructVectorDestroy,(ex->ss_b));
1100   PetscFunctionReturn(0);
1101 }
1102 
1103 EXTERN_C_BEGIN
1104 #undef __FUNCT__
1105 #define __FUNCT__ "MatCreate_HYPRESStruct"
1106 PetscErrorCode  MatCreate_HYPRESStruct(Mat B)
1107 {
1108   Mat_HYPRESStruct *ex;
1109   PetscErrorCode   ierr;
1110 
1111   PetscFunctionBegin;
1112   ierr            = PetscNewLog(B,Mat_HYPRESStruct,&ex);CHKERRQ(ierr);
1113   B->data         = (void*)ex;
1114   B->rmap->bs     = 1;
1115   B->assembled    = PETSC_FALSE;
1116 
1117   B->insertmode   = NOT_SET_VALUES;
1118 
1119   B->ops->assemblyend    = MatAssemblyEnd_HYPRESStruct;
1120   B->ops->mult           = MatMult_HYPRESStruct;
1121   B->ops->zeroentries    = MatZeroEntries_HYPRESStruct;
1122   B->ops->destroy        = MatDestroy_HYPRESStruct;
1123 
1124   ex->needsinitialization = PETSC_TRUE;
1125 
1126   ierr = MPI_Comm_dup(((PetscObject)B)->comm,&(ex->hcomm));CHKERRQ(ierr);
1127   ierr = PetscObjectComposeFunctionDynamic((PetscObject)B,"MatSetDM_C","MatSetDM_HYPRESStruct",MatSetDM_HYPRESStruct);CHKERRQ(ierr);
1128   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRESSTRUCT);CHKERRQ(ierr);
1129   PetscFunctionReturn(0);
1130 }
1131 EXTERN_C_END
1132 
1133 
1134 
1135