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