xref: /petsc/src/dm/impls/da/hypre/mhyp.c (revision f6dae198d6c3cf07add646cdb0a134c2f8c0b9f9)
1f4bdf6c4SBarry Smith 
2f4bdf6c4SBarry Smith /*
3f4bdf6c4SBarry Smith     Creates hypre ijmatrix from PETSc matrix
4f4bdf6c4SBarry Smith */
5c6db04a5SJed Brown #include <petscsys.h>
6af0996ceSBarry Smith #include <petsc/private/matimpl.h>
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;
151a83f524SJed Brown   PetscInt       i,n_d,n_o;
161a83f524SJed Brown   const PetscInt *ia_d,*ia_o;
17f4bdf6c4SBarry Smith   PetscBool      done_d=PETSC_FALSE,done_o=PETSC_FALSE;
180298fd71SBarry Smith   PetscInt       *nnz_d=NULL,*nnz_o=NULL;
19f4bdf6c4SBarry Smith 
20f4bdf6c4SBarry Smith   PetscFunctionBegin;
21f4bdf6c4SBarry Smith   if (A_d) { /* determine number of nonzero entries in local diagonal part */
220298fd71SBarry Smith     ierr = MatGetRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,&n_d,&ia_d,NULL,&done_d);CHKERRQ(ierr);
23f4bdf6c4SBarry Smith     if (done_d) {
24785e854fSJed Brown       ierr = PetscMalloc1(n_d,&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     }
299383aa05SJed Brown     ierr = MatRestoreRowIJ(A_d,0,PETSC_FALSE,PETSC_FALSE,NULL,&ia_d,NULL,&done_d);CHKERRQ(ierr);
30f4bdf6c4SBarry Smith   }
31f4bdf6c4SBarry Smith   if (A_o) { /* determine number of nonzero entries in local off-diagonal part */
320298fd71SBarry Smith     ierr = MatGetRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,NULL,&done_o);CHKERRQ(ierr);
33f4bdf6c4SBarry Smith     if (done_o) {
34785e854fSJed Brown       ierr = PetscMalloc1(n_o,&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     }
390298fd71SBarry Smith     ierr = MatRestoreRowIJ(A_o,0,PETSC_FALSE,PETSC_FALSE,&n_o,&ia_o,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 */
43785e854fSJed Brown       ierr = PetscMalloc1(n_d,&nnz_o);CHKERRQ(ierr);
44f4bdf6c4SBarry Smith       for (i=0; i<n_d; i++) {
45f4bdf6c4SBarry Smith         nnz_o[i] = 0;
46f4bdf6c4SBarry Smith       }
47f4bdf6c4SBarry Smith     }
488b1f7689SBarry Smith     PetscStackCallStandard(HYPRE_IJMatrixSetDiagOffdSizes,(ij,(HYPRE_Int *)nnz_d,(HYPRE_Int *)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;
604ddd07fcSJed Brown   PetscInt       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);
664994cf47SJed Brown   ierr   = MatSetUp(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;
71ce94432eSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij));
72fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
73f4bdf6c4SBarry Smith   {
74f4bdf6c4SBarry Smith     PetscBool      same;
75f4bdf6c4SBarry Smith     Mat            A_d,A_o;
761a83f524SJed Brown     const PetscInt *colmap;
77251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((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     }
83251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((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     }
89251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&same);CHKERRQ(ierr);
90f4bdf6c4SBarry Smith     if (same) {
910298fd71SBarry Smith       ierr = MatHYPRE_IJMatrixPreallocate(A,NULL,*ij);CHKERRQ(ierr);
92f4bdf6c4SBarry Smith       PetscFunctionReturn(0);
93f4bdf6c4SBarry Smith     }
94251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&same);CHKERRQ(ierr);
95f4bdf6c4SBarry Smith     if (same) {
960298fd71SBarry Smith       ierr = MatHYPRE_IJMatrixPreallocate(A,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;
1144cb006feSStefano Zampini   PetscInt          i,rstart,rend,ncols,nr,nc;
115f4bdf6c4SBarry Smith   const PetscScalar *values;
116f4bdf6c4SBarry Smith   const PetscInt    *cols;
117f4bdf6c4SBarry Smith   PetscBool         flg;
118f4bdf6c4SBarry Smith 
119f4bdf6c4SBarry Smith   PetscFunctionBegin;
120251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
1214cb006feSStefano Zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
1224cb006feSStefano Zampini   if (flg && nr == nc) {
123f4bdf6c4SBarry Smith     ierr = MatHYPRE_IJMatrixFastCopy_MPIAIJ(A,ij);CHKERRQ(ierr);
124f4bdf6c4SBarry Smith     PetscFunctionReturn(0);
125f4bdf6c4SBarry Smith   }
126251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQAIJ,&flg);CHKERRQ(ierr);
127f4bdf6c4SBarry Smith   if (flg) {
128f4bdf6c4SBarry Smith     ierr = MatHYPRE_IJMatrixFastCopy_SeqAIJ(A,ij);CHKERRQ(ierr);
129f4bdf6c4SBarry Smith     PetscFunctionReturn(0);
130f4bdf6c4SBarry Smith   }
131f4bdf6c4SBarry Smith 
132f4bdf6c4SBarry Smith   ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
133fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
134f4bdf6c4SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
135f4bdf6c4SBarry Smith   for (i=rstart; i<rend; i++) {
136f4bdf6c4SBarry Smith     ierr = MatGetRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
1378b1f7689SBarry Smith     PetscStackCallStandard(HYPRE_IJMatrixSetValues,(ij,1,(HYPRE_Int *)&ncols,(HYPRE_Int *)&i,(HYPRE_Int *)cols,values));
138f4bdf6c4SBarry Smith     ierr = MatRestoreRow(A,i,&ncols,&cols,&values);CHKERRQ(ierr);
139f4bdf6c4SBarry Smith   }
140fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(ij));
141f4bdf6c4SBarry Smith   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
142f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
143f4bdf6c4SBarry Smith }
144f4bdf6c4SBarry Smith 
145f4bdf6c4SBarry Smith /*
146f4bdf6c4SBarry Smith     This copies the CSR format directly from the PETSc data structure to the hypre
147f4bdf6c4SBarry Smith     data structure without calls to MatGetRow() or hypre's set values.
148f4bdf6c4SBarry Smith 
149f4bdf6c4SBarry Smith */
150c6db04a5SJed Brown #include <_hypre_IJ_mv.h>
151c6db04a5SJed Brown #include <HYPRE_IJ_mv.h>
152c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
153f4bdf6c4SBarry Smith 
154f4bdf6c4SBarry Smith #undef __FUNCT__
155f4bdf6c4SBarry Smith #define __FUNCT__ "MatHYPRE_IJMatrixFastCopy_SeqAIJ"
156f4bdf6c4SBarry Smith PetscErrorCode MatHYPRE_IJMatrixFastCopy_SeqAIJ(Mat A,HYPRE_IJMatrix ij)
157f4bdf6c4SBarry Smith {
158f4bdf6c4SBarry Smith   PetscErrorCode ierr;
159e497e08cSBarry Smith   Mat_SeqAIJ     *pdiag = (Mat_SeqAIJ*)A->data;
160f4bdf6c4SBarry Smith 
161f4bdf6c4SBarry Smith   hypre_ParCSRMatrix    *par_matrix;
162f4bdf6c4SBarry Smith   hypre_AuxParCSRMatrix *aux_matrix;
16325578ef6SJed Brown   hypre_CSRMatrix       *hdiag;
164f4bdf6c4SBarry Smith 
165f4bdf6c4SBarry Smith   PetscFunctionBegin;
166f4bdf6c4SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
167f4bdf6c4SBarry Smith   PetscValidType(A,1);
168f4bdf6c4SBarry Smith   PetscValidPointer(ij,2);
169f4bdf6c4SBarry Smith 
170f4bdf6c4SBarry Smith   ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
171fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(ij));
172f4bdf6c4SBarry Smith   par_matrix = (hypre_ParCSRMatrix*)hypre_IJMatrixObject(ij);
173f4bdf6c4SBarry Smith   aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(ij);
174f4bdf6c4SBarry Smith   hdiag      = hypre_ParCSRMatrixDiag(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;
184fd3f9acdSBarry Smith   PetscStackCallStandard(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 */
2090298fd71SBarry Smith   ierr = MatGetOwnershipRange(A,&cstart,NULL);CHKERRQ(ierr);
210f4bdf6c4SBarry Smith 
211f4bdf6c4SBarry Smith   ierr = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
212fd3f9acdSBarry Smith   PetscStackCallStandard(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 */
2238b1f7689SBarry Smith   jj  = (PetscInt*)hdiag->j;
2248b1f7689SBarry Smith   pjj = (PetscInt*)pdiag->j;
2258865f1eaSKarl Rupp   for (i=0; i<pdiag->nz; i++) jj[i] = cstart + pjj[i];
226f4bdf6c4SBarry Smith   ierr = PetscMemcpy(hdiag->data,pdiag->a,pdiag->nz*sizeof(PetscScalar));
227f4bdf6c4SBarry Smith   ierr = PetscMemcpy(hoffd->i,poffd->i,(pA->A->rmap->n + 1)*sizeof(PetscInt));
228f4bdf6c4SBarry Smith   /* need to move the offd column indices (hoffd->j) back to global numbering since hypre is expecting this
229f4bdf6c4SBarry Smith      If we hacked a hypre a bit more we might be able to avoid this step */
230b0de412dSBarry Smith   jj  = (PetscInt*) hoffd->j;
231b0de412dSBarry Smith   pjj = (PetscInt*) poffd->j;
2328865f1eaSKarl Rupp   for (i=0; i<poffd->nz; i++) jj[i] = garray[pjj[i]];
233f4bdf6c4SBarry Smith   ierr = PetscMemcpy(hoffd->data,poffd->a,poffd->nz*sizeof(PetscScalar));
234f4bdf6c4SBarry Smith 
235f4bdf6c4SBarry Smith   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
236fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(ij));
237f4bdf6c4SBarry Smith   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
238f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
239f4bdf6c4SBarry Smith }
240f4bdf6c4SBarry Smith 
241f4bdf6c4SBarry Smith /*
242f4bdf6c4SBarry Smith     Does NOT copy the data over, instead uses DIRECTLY the pointers from the PETSc MPIAIJ format
243f4bdf6c4SBarry Smith 
244f4bdf6c4SBarry Smith     This is UNFINISHED and does NOT work! The problem is that hypre puts the diagonal entry first
245f4bdf6c4SBarry Smith     which will corrupt the PETSc data structure if we did this. Need a work around to this problem.
246f4bdf6c4SBarry Smith */
247c6db04a5SJed Brown #include <_hypre_IJ_mv.h>
248c6db04a5SJed Brown #include <HYPRE_IJ_mv.h>
249f4bdf6c4SBarry Smith 
250f4bdf6c4SBarry Smith #undef __FUNCT__
251f4bdf6c4SBarry Smith #define __FUNCT__ "MatHYPRE_IJMatrixLink"
252f4bdf6c4SBarry Smith PetscErrorCode MatHYPRE_IJMatrixLink(Mat A,HYPRE_IJMatrix *ij)
253f4bdf6c4SBarry Smith {
254f4bdf6c4SBarry Smith   PetscErrorCode        ierr;
2554ddd07fcSJed Brown   PetscInt              rstart,rend,cstart,cend;
256f4bdf6c4SBarry Smith   PetscBool             flg;
257f4bdf6c4SBarry Smith   hypre_AuxParCSRMatrix *aux_matrix;
258f4bdf6c4SBarry Smith 
259f4bdf6c4SBarry Smith   PetscFunctionBegin;
260f4bdf6c4SBarry Smith   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
261f4bdf6c4SBarry Smith   PetscValidType(A,1);
262f4bdf6c4SBarry Smith   PetscValidPointer(ij,2);
263251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);CHKERRQ(ierr);
264ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"Can only use with PETSc MPIAIJ matrices");
2654994cf47SJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
266f4bdf6c4SBarry Smith 
267f4bdf6c4SBarry Smith   ierr   = PetscLogEventBegin(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
268f4bdf6c4SBarry Smith   rstart = A->rmap->rstart;
269f4bdf6c4SBarry Smith   rend   = A->rmap->rend;
270f4bdf6c4SBarry Smith   cstart = A->cmap->rstart;
271f4bdf6c4SBarry Smith   cend   = A->cmap->rend;
272ce94432eSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixCreate,(PetscObjectComm((PetscObject)A),rstart,rend-1,cstart,cend-1,ij));
273fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixSetObjectType,(*ij,HYPRE_PARCSR));
274f4bdf6c4SBarry Smith 
275fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixInitialize,(*ij));
276fd3f9acdSBarry Smith   PetscStackCall("hypre_IJMatrixTranslator",aux_matrix = (hypre_AuxParCSRMatrix*)hypre_IJMatrixTranslator(*ij));
277f4bdf6c4SBarry Smith 
278f4bdf6c4SBarry Smith   hypre_AuxParCSRMatrixNeedAux(aux_matrix) = 0;
279f4bdf6c4SBarry Smith 
280f4bdf6c4SBarry Smith   /* this is the Hack part where we monkey directly with the hypre datastructures */
281f4bdf6c4SBarry Smith 
282fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_IJMatrixAssemble,(*ij));
283f4bdf6c4SBarry Smith   ierr = PetscLogEventEnd(MAT_Convert,A,0,0,0);CHKERRQ(ierr);
284f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
285f4bdf6c4SBarry Smith }
286f4bdf6c4SBarry Smith 
287f4bdf6c4SBarry Smith /* -----------------------------------------------------------------------------------------------------------------*/
288f4bdf6c4SBarry Smith 
289f4bdf6c4SBarry Smith /*MC
290f4bdf6c4SBarry Smith    MATHYPRESTRUCT - MATHYPRESTRUCT = "hyprestruct" - A matrix type to be used for parallel sparse matrices
291f4bdf6c4SBarry Smith           based on the hypre HYPRE_StructMatrix.
292f4bdf6c4SBarry Smith 
293f4bdf6c4SBarry Smith    Level: intermediate
294f4bdf6c4SBarry Smith 
295f4bdf6c4SBarry Smith    Notes: Unlike the more general support for blocks in hypre this allows only one block per process and requires the block
296f4bdf6c4SBarry Smith           be defined by a DMDA.
297f4bdf6c4SBarry Smith 
298c688c046SMatthew G Knepley           The matrix needs a DMDA associated with it by either a call to MatSetupDM() or if the matrix is obtained from DMCreateMatrix()
299f4bdf6c4SBarry Smith 
300c688c046SMatthew G Knepley .seealso: MatCreate(), PCPFMG, MatSetupDM(), DMCreateMatrix()
301f4bdf6c4SBarry Smith M*/
302f4bdf6c4SBarry Smith 
303f4bdf6c4SBarry Smith 
304f4bdf6c4SBarry Smith #undef __FUNCT__
305f4bdf6c4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_HYPREStruct_3d"
306f4bdf6c4SBarry Smith PetscErrorCode  MatSetValuesLocal_HYPREStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
307f4bdf6c4SBarry Smith {
308f4bdf6c4SBarry Smith   PetscErrorCode    ierr;
309f4bdf6c4SBarry Smith   PetscInt          i,j,stencil,index[3],row,entries[7];
310f4bdf6c4SBarry Smith   const PetscScalar *values = y;
311f4bdf6c4SBarry Smith   Mat_HYPREStruct   *ex     = (Mat_HYPREStruct*) mat->data;
312f4bdf6c4SBarry Smith 
313f4bdf6c4SBarry Smith   PetscFunctionBegin;
314*f6dae198SJed Brown   if (PetscUnlikely(ncol >= 7)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ncol %D >= 7 too large",ncol);
315f4bdf6c4SBarry Smith   for (i=0; i<nrow; i++) {
316f4bdf6c4SBarry Smith     for (j=0; j<ncol; j++) {
317f4bdf6c4SBarry Smith       stencil = icol[j] - irow[i];
318f4bdf6c4SBarry Smith       if (!stencil) {
319f4bdf6c4SBarry Smith         entries[j] = 3;
320f4bdf6c4SBarry Smith       } else if (stencil == -1) {
321f4bdf6c4SBarry Smith         entries[j] = 2;
322f4bdf6c4SBarry Smith       } else if (stencil == 1) {
323f4bdf6c4SBarry Smith         entries[j] = 4;
324f4bdf6c4SBarry Smith       } else if (stencil == -ex->gnx) {
325f4bdf6c4SBarry Smith         entries[j] = 1;
326f4bdf6c4SBarry Smith       } else if (stencil == ex->gnx) {
327f4bdf6c4SBarry Smith         entries[j] = 5;
328f4bdf6c4SBarry Smith       } else if (stencil == -ex->gnxgny) {
329f4bdf6c4SBarry Smith         entries[j] = 0;
330f4bdf6c4SBarry Smith       } else if (stencil == ex->gnxgny) {
331f4bdf6c4SBarry Smith         entries[j] = 6;
332f4bdf6c4SBarry 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);
333f4bdf6c4SBarry Smith     }
334f4bdf6c4SBarry Smith     row      = ex->gindices[irow[i]] - ex->rstart;
335f4bdf6c4SBarry Smith     index[0] = ex->xs + (row % ex->nx);
336f4bdf6c4SBarry Smith     index[1] = ex->ys + ((row/ex->nx) % ex->ny);
337f4bdf6c4SBarry Smith     index[2] = ex->zs + (row/(ex->nxny));
338f4bdf6c4SBarry Smith     if (addv == ADD_VALUES) {
3398b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_StructMatrixAddToValues,(ex->hmat,(HYPRE_Int *)index,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
340f4bdf6c4SBarry Smith     } else {
3418b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_StructMatrixSetValues,(ex->hmat,(HYPRE_Int *)index,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
342f4bdf6c4SBarry Smith     }
343f4bdf6c4SBarry Smith     values += ncol;
344f4bdf6c4SBarry Smith   }
345f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
346f4bdf6c4SBarry Smith }
347f4bdf6c4SBarry Smith 
348f4bdf6c4SBarry Smith #undef __FUNCT__
349f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_HYPREStruct_3d"
350f4bdf6c4SBarry Smith PetscErrorCode  MatZeroRowsLocal_HYPREStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscScalar d,Vec x,Vec b)
351f4bdf6c4SBarry Smith {
352f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
353f4bdf6c4SBarry Smith   PetscInt        i,index[3],row,entries[7] = {0,1,2,3,4,5,6};
354f4bdf6c4SBarry Smith   PetscScalar     values[7];
355f4bdf6c4SBarry Smith   Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data;
356f4bdf6c4SBarry Smith 
357f4bdf6c4SBarry Smith   PetscFunctionBegin;
358ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No support");
359f4bdf6c4SBarry Smith   ierr      = PetscMemzero(values,7*sizeof(PetscScalar));CHKERRQ(ierr);
360f4bdf6c4SBarry Smith   values[3] = d;
361f4bdf6c4SBarry Smith   for (i=0; i<nrow; i++) {
362f4bdf6c4SBarry Smith     row      = ex->gindices[irow[i]] - ex->rstart;
363f4bdf6c4SBarry Smith     index[0] = ex->xs + (row % ex->nx);
364f4bdf6c4SBarry Smith     index[1] = ex->ys + ((row/ex->nx) % ex->ny);
365f4bdf6c4SBarry Smith     index[2] = ex->zs + (row/(ex->nxny));
3668b1f7689SBarry Smith     PetscStackCallStandard(HYPRE_StructMatrixSetValues,(ex->hmat,(HYPRE_Int *)index,7,(HYPRE_Int *)entries,values));
367f4bdf6c4SBarry Smith   }
368fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixAssemble,(ex->hmat));
369f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
370f4bdf6c4SBarry Smith }
371f4bdf6c4SBarry Smith 
372f4bdf6c4SBarry Smith #undef __FUNCT__
373f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPREStruct_3d"
374f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPREStruct_3d(Mat mat)
375f4bdf6c4SBarry Smith {
376f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
377f4bdf6c4SBarry Smith   PetscInt        indices[7] = {0,1,2,3,4,5,6};
378f4bdf6c4SBarry Smith   Mat_HYPREStruct *ex        = (Mat_HYPREStruct*) mat->data;
379f4bdf6c4SBarry Smith 
380f4bdf6c4SBarry Smith   PetscFunctionBegin;
381f4bdf6c4SBarry Smith   /* hypre has no public interface to do this */
3828b1f7689SBarry Smith   PetscStackCallStandard(hypre_StructMatrixClearBoxValues,(ex->hmat,&ex->hbox,7,(HYPRE_Int *)indices,0,1));
383fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixAssemble,(ex->hmat));
384f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
385f4bdf6c4SBarry Smith }
386f4bdf6c4SBarry Smith 
387f4bdf6c4SBarry Smith #undef __FUNCT__
388c688c046SMatthew G Knepley #define __FUNCT__ "MatSetupDM_HYPREStruct"
3891e6b0712SBarry Smith static PetscErrorCode  MatSetupDM_HYPREStruct(Mat mat,DM da)
390f4bdf6c4SBarry Smith {
391f4bdf6c4SBarry Smith   PetscErrorCode         ierr;
392f4bdf6c4SBarry Smith   Mat_HYPREStruct        *ex = (Mat_HYPREStruct*) mat->data;
3934ddd07fcSJed Brown   PetscInt               dim,dof,sw[3],nx,ny,nz,ilower[3],iupper[3],ssize,i;
394bff4a2f0SMatthew G. Knepley   DMBoundaryType         px,py,pz;
395f4bdf6c4SBarry Smith   DMDAStencilType        st;
396565245c5SBarry Smith   ISLocalToGlobalMapping ltog;
397f4bdf6c4SBarry Smith 
398f4bdf6c4SBarry Smith   PetscFunctionBegin;
399f4bdf6c4SBarry Smith   ex->da = da;
400f4bdf6c4SBarry Smith   ierr   = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);
401f4bdf6c4SBarry Smith 
4027ab44359SJed Brown   ierr       = DMDAGetInfo(ex->da,&dim,0,0,0,0,0,0,&dof,&sw[0],&px,&py,&pz,&st);CHKERRQ(ierr);
403f4bdf6c4SBarry Smith   ierr       = DMDAGetCorners(ex->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
404f4bdf6c4SBarry Smith   iupper[0] += ilower[0] - 1;
405f4bdf6c4SBarry Smith   iupper[1] += ilower[1] - 1;
406f4bdf6c4SBarry Smith   iupper[2] += ilower[2] - 1;
407f4bdf6c4SBarry Smith 
408f4bdf6c4SBarry Smith   /* the hypre_Box is used to zero out the matrix entries in MatZeroValues() */
409f4bdf6c4SBarry Smith   ex->hbox.imin[0] = ilower[0];
410f4bdf6c4SBarry Smith   ex->hbox.imin[1] = ilower[1];
411f4bdf6c4SBarry Smith   ex->hbox.imin[2] = ilower[2];
412f4bdf6c4SBarry Smith   ex->hbox.imax[0] = iupper[0];
413f4bdf6c4SBarry Smith   ex->hbox.imax[1] = iupper[1];
414f4bdf6c4SBarry Smith   ex->hbox.imax[2] = iupper[2];
415f4bdf6c4SBarry Smith 
416f4bdf6c4SBarry Smith   /* create the hypre grid object and set its information */
417ce94432eSBarry Smith   if (dof > 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Currently only support for scalar problems");
418ce94432eSBarry Smith   if (px || py || pz) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add periodic support by calling HYPRE_StructGridSetPeriodic()");
419fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructGridCreate,(ex->hcomm,dim,&ex->hgrid));
4208b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_StructGridSetExtents,(ex->hgrid,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper));
421fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructGridAssemble,(ex->hgrid));
422f4bdf6c4SBarry Smith 
423f4bdf6c4SBarry Smith   sw[1] = sw[0];
424f4bdf6c4SBarry Smith   sw[2] = sw[1];
4258b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_StructGridSetNumGhost,(ex->hgrid,(HYPRE_Int *)sw));
426f4bdf6c4SBarry Smith 
427f4bdf6c4SBarry Smith   /* create the hypre stencil object and set its information */
428ce94432eSBarry Smith   if (sw[0] > 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add support for wider stencils");
429ce94432eSBarry Smith   if (st == DMDA_STENCIL_BOX) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add support for box stencils");
430f4bdf6c4SBarry Smith   if (dim == 1) {
4314ddd07fcSJed Brown     PetscInt offsets[3][1] = {{-1},{0},{1}};
432f4bdf6c4SBarry Smith     ssize = 3;
433fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil));
434f4bdf6c4SBarry Smith     for (i=0; i<ssize; i++) {
4358b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_StructStencilSetElement,(ex->hstencil,i,(HYPRE_Int *)offsets[i]));
436f4bdf6c4SBarry Smith     }
437f4bdf6c4SBarry Smith   } else if (dim == 2) {
4384ddd07fcSJed Brown     PetscInt offsets[5][2] = {{0,-1},{-1,0},{0,0},{1,0},{0,1}};
439f4bdf6c4SBarry Smith     ssize = 5;
440fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil));
441f4bdf6c4SBarry Smith     for (i=0; i<ssize; i++) {
4428b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_StructStencilSetElement,(ex->hstencil,i,(HYPRE_Int *)offsets[i]));
443f4bdf6c4SBarry Smith     }
444f4bdf6c4SBarry Smith   } else if (dim == 3) {
4454ddd07fcSJed Brown     PetscInt offsets[7][3] = {{0,0,-1},{0,-1,0},{-1,0,0},{0,0,0},{1,0,0},{0,1,0},{0,0,1}};
446f4bdf6c4SBarry Smith     ssize = 7;
447fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil));
448f4bdf6c4SBarry Smith     for (i=0; i<ssize; i++) {
4498b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_StructStencilSetElement,(ex->hstencil,i,(HYPRE_Int *)offsets[i]));
450f4bdf6c4SBarry Smith     }
451f4bdf6c4SBarry Smith   }
452f4bdf6c4SBarry Smith 
453f4bdf6c4SBarry Smith   /* create the HYPRE vector for rhs and solution */
454fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorCreate,(ex->hcomm,ex->hgrid,&ex->hb));
455fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorCreate,(ex->hcomm,ex->hgrid,&ex->hx));
456fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorInitialize,(ex->hb));
457fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorInitialize,(ex->hx));
458fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorAssemble,(ex->hb));
459fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorAssemble,(ex->hx));
460f4bdf6c4SBarry Smith 
461f4bdf6c4SBarry Smith   /* create the hypre matrix object and set its information */
462fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixCreate,(ex->hcomm,ex->hgrid,ex->hstencil,&ex->hmat));
463fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructGridDestroy,(ex->hgrid));
464fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructStencilDestroy,(ex->hstencil));
465f4bdf6c4SBarry Smith   if (ex->needsinitialization) {
466fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructMatrixInitialize,(ex->hmat));
467f4bdf6c4SBarry Smith     ex->needsinitialization = PETSC_FALSE;
468f4bdf6c4SBarry Smith   }
469f4bdf6c4SBarry Smith 
470f4bdf6c4SBarry Smith   /* set the global and local sizes of the matrix */
471f4bdf6c4SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
472f4bdf6c4SBarry Smith   ierr = MatSetSizes(mat,dof*nx*ny*nz,dof*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
473f4bdf6c4SBarry Smith   ierr = PetscLayoutSetBlockSize(mat->rmap,1);CHKERRQ(ierr);
474f4bdf6c4SBarry Smith   ierr = PetscLayoutSetBlockSize(mat->cmap,1);CHKERRQ(ierr);
475f4bdf6c4SBarry Smith   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
476f4bdf6c4SBarry Smith   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
477f4bdf6c4SBarry Smith 
478f4bdf6c4SBarry Smith   if (dim == 3) {
479f4bdf6c4SBarry Smith     mat->ops->setvalueslocal = MatSetValuesLocal_HYPREStruct_3d;
480f4bdf6c4SBarry Smith     mat->ops->zerorowslocal  = MatZeroRowsLocal_HYPREStruct_3d;
481f4bdf6c4SBarry Smith     mat->ops->zeroentries    = MatZeroEntries_HYPREStruct_3d;
4828865f1eaSKarl Rupp 
483f4bdf6c4SBarry Smith     ierr = MatZeroEntries_HYPREStruct_3d(mat);CHKERRQ(ierr);
484ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only support for 3d DMDA currently");
485f4bdf6c4SBarry Smith 
486f4bdf6c4SBarry Smith   /* get values that will be used repeatedly in MatSetValuesLocal() and MatZeroRowsLocal() repeatedly */
4870298fd71SBarry Smith   ierr        = MatGetOwnershipRange(mat,&ex->rstart,NULL);CHKERRQ(ierr);
488565245c5SBarry Smith   ierr        = DMGetLocalToGlobalMapping(ex->da,&ltog);CHKERRQ(ierr);
489565245c5SBarry Smith   ierr        = ISLocalToGlobalMappingGetIndices(ltog, (const PetscInt **) &ex->gindices);CHKERRQ(ierr);
490f4bdf6c4SBarry Smith   ierr        = DMDAGetGhostCorners(ex->da,0,0,0,&ex->gnx,&ex->gnxgny,0);CHKERRQ(ierr);
491f4bdf6c4SBarry Smith   ex->gnxgny *= ex->gnx;
492f4bdf6c4SBarry Smith   ierr        = DMDAGetCorners(ex->da,&ex->xs,&ex->ys,&ex->zs,&ex->nx,&ex->ny,0);CHKERRQ(ierr);
493f4bdf6c4SBarry Smith   ex->nxny    = ex->nx*ex->ny;
494f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
495f4bdf6c4SBarry Smith }
496f4bdf6c4SBarry Smith 
497f4bdf6c4SBarry Smith #undef __FUNCT__
498f4bdf6c4SBarry Smith #define __FUNCT__ "MatMult_HYPREStruct"
499f4bdf6c4SBarry Smith PetscErrorCode MatMult_HYPREStruct(Mat A,Vec x,Vec y)
500f4bdf6c4SBarry Smith {
501f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
502f4bdf6c4SBarry Smith   PetscScalar     *xx,*yy;
5034ddd07fcSJed Brown   PetscInt        ilower[3],iupper[3];
504f4bdf6c4SBarry Smith   Mat_HYPREStruct *mx = (Mat_HYPREStruct*)(A->data);
505f4bdf6c4SBarry Smith 
506f4bdf6c4SBarry Smith   PetscFunctionBegin;
507f4bdf6c4SBarry Smith   ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
5088865f1eaSKarl Rupp 
509f4bdf6c4SBarry Smith   iupper[0] += ilower[0] - 1;
510f4bdf6c4SBarry Smith   iupper[1] += ilower[1] - 1;
511f4bdf6c4SBarry Smith   iupper[2] += ilower[2] - 1;
512f4bdf6c4SBarry Smith 
513f4bdf6c4SBarry Smith   /* copy x values over to hypre */
514fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
515f4bdf6c4SBarry Smith   ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
5168b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_StructVectorSetBoxValues,(mx->hb,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,xx));
517f4bdf6c4SBarry Smith   ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
518fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorAssemble,(mx->hb));
519fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixMatvec,(1.0,mx->hmat,mx->hb,0.0,mx->hx));
520f4bdf6c4SBarry Smith 
521f4bdf6c4SBarry Smith   /* copy solution values back to PETSc */
522f4bdf6c4SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
5238b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_StructVectorGetBoxValues,(mx->hx,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,yy));
524f4bdf6c4SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
525f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
526f4bdf6c4SBarry Smith }
527f4bdf6c4SBarry Smith 
528f4bdf6c4SBarry Smith #undef __FUNCT__
529f4bdf6c4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_HYPREStruct"
530f4bdf6c4SBarry Smith PetscErrorCode MatAssemblyEnd_HYPREStruct(Mat mat,MatAssemblyType mode)
531f4bdf6c4SBarry Smith {
532f4bdf6c4SBarry Smith   Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data;
533f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
534f4bdf6c4SBarry Smith 
535f4bdf6c4SBarry Smith   PetscFunctionBegin;
536fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixAssemble,(ex->hmat));
537fd3f9acdSBarry Smith   /* PetscStackCallStandard(HYPRE_StructMatrixPrint,("dummy",ex->hmat,0)); */
538f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
539f4bdf6c4SBarry Smith }
540f4bdf6c4SBarry Smith 
541f4bdf6c4SBarry Smith #undef __FUNCT__
542f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPREStruct"
543f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPREStruct(Mat mat)
544f4bdf6c4SBarry Smith {
545f4bdf6c4SBarry Smith   PetscFunctionBegin;
546f4bdf6c4SBarry Smith   /* before the DMDA is set to the matrix the zero doesn't need to do anything */
547f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
548f4bdf6c4SBarry Smith }
549f4bdf6c4SBarry Smith 
550f4bdf6c4SBarry Smith 
551f4bdf6c4SBarry Smith #undef __FUNCT__
552f4bdf6c4SBarry Smith #define __FUNCT__ "MatDestroy_HYPREStruct"
553f4bdf6c4SBarry Smith PetscErrorCode MatDestroy_HYPREStruct(Mat mat)
554f4bdf6c4SBarry Smith {
555f4bdf6c4SBarry Smith   Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data;
556f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
557f4bdf6c4SBarry Smith 
558f4bdf6c4SBarry Smith   PetscFunctionBegin;
559fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixDestroy,(ex->hmat));
560fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorDestroy,(ex->hx));
561fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorDestroy,(ex->hb));
562f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
563f4bdf6c4SBarry Smith }
564f4bdf6c4SBarry Smith 
565f4bdf6c4SBarry Smith 
566f4bdf6c4SBarry Smith #undef __FUNCT__
567f4bdf6c4SBarry Smith #define __FUNCT__ "MatCreate_HYPREStruct"
5688cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_HYPREStruct(Mat B)
569f4bdf6c4SBarry Smith {
570f4bdf6c4SBarry Smith   Mat_HYPREStruct *ex;
571f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
572f4bdf6c4SBarry Smith 
573f4bdf6c4SBarry Smith   PetscFunctionBegin;
574b00a9115SJed Brown   ierr         = PetscNewLog(B,&ex);CHKERRQ(ierr);
575f4bdf6c4SBarry Smith   B->data      = (void*)ex;
576f4bdf6c4SBarry Smith   B->rmap->bs  = 1;
577f4bdf6c4SBarry Smith   B->assembled = PETSC_FALSE;
578f4bdf6c4SBarry Smith 
579f4bdf6c4SBarry Smith   B->insertmode = NOT_SET_VALUES;
580f4bdf6c4SBarry Smith 
581f4bdf6c4SBarry Smith   B->ops->assemblyend = MatAssemblyEnd_HYPREStruct;
582f4bdf6c4SBarry Smith   B->ops->mult        = MatMult_HYPREStruct;
583f4bdf6c4SBarry Smith   B->ops->zeroentries = MatZeroEntries_HYPREStruct;
584f4bdf6c4SBarry Smith   B->ops->destroy     = MatDestroy_HYPREStruct;
585f4bdf6c4SBarry Smith 
586f4bdf6c4SBarry Smith   ex->needsinitialization = PETSC_TRUE;
587f4bdf6c4SBarry Smith 
588ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&(ex->hcomm));CHKERRQ(ierr);
589bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSetupDM_C",MatSetupDM_HYPREStruct);CHKERRQ(ierr);
590f4bdf6c4SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRESTRUCT);CHKERRQ(ierr);
591f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
592f4bdf6c4SBarry Smith }
593f4bdf6c4SBarry Smith 
594f4bdf6c4SBarry Smith /*MC
595f4bdf6c4SBarry Smith    MATHYPRESSTRUCT - MATHYPRESSTRUCT = "hypresstruct" - A matrix type to be used for parallel sparse matrices
596f4bdf6c4SBarry Smith           based on the hypre HYPRE_SStructMatrix.
597f4bdf6c4SBarry Smith 
598f4bdf6c4SBarry Smith 
599f4bdf6c4SBarry Smith    Level: intermediate
600f4bdf6c4SBarry Smith 
601f4bdf6c4SBarry Smith    Notes: Unlike hypre's general semi-struct object consisting of a collection of structured-grid objects and unstructured
602f4bdf6c4SBarry Smith           grid objects, we will restrict the semi-struct objects to consist of only structured-grid components.
603f4bdf6c4SBarry Smith 
604f4bdf6c4SBarry 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
605f4bdf6c4SBarry Smith           be defined by a DMDA.
606f4bdf6c4SBarry Smith 
607c688c046SMatthew G Knepley           The matrix needs a DMDA associated with it by either a call to MatSetupDM() or if the matrix is obtained from DMCreateMatrix()
608f4bdf6c4SBarry Smith 
609f4bdf6c4SBarry Smith M*/
610f4bdf6c4SBarry Smith 
611f4bdf6c4SBarry Smith #undef __FUNCT__
612f4bdf6c4SBarry Smith #define __FUNCT__ "MatSetValuesLocal_HYPRESStruct_3d"
613f4bdf6c4SBarry Smith PetscErrorCode  MatSetValuesLocal_HYPRESStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
614f4bdf6c4SBarry Smith {
615f4bdf6c4SBarry Smith   PetscErrorCode    ierr;
616f4bdf6c4SBarry Smith   PetscInt          i,j,stencil,index[3];
617f4bdf6c4SBarry Smith   const PetscScalar *values = y;
618f4bdf6c4SBarry Smith   Mat_HYPRESStruct  *ex     = (Mat_HYPRESStruct*) mat->data;
619f4bdf6c4SBarry Smith 
6204ddd07fcSJed Brown   PetscInt part = 0;          /* Petsc sstruct interface only allows 1 part */
6214ddd07fcSJed Brown   PetscInt ordering;
6224ddd07fcSJed Brown   PetscInt grid_rank, to_grid_rank;
6234ddd07fcSJed Brown   PetscInt var_type, to_var_type;
6244ddd07fcSJed Brown   PetscInt to_var_entry = 0;
625f4bdf6c4SBarry Smith 
6264ddd07fcSJed Brown   PetscInt nvars= ex->nvars;
627f4bdf6c4SBarry Smith   PetscInt row,*entries;
628f4bdf6c4SBarry Smith 
629f4bdf6c4SBarry Smith   PetscFunctionBegin;
630785e854fSJed Brown   ierr    = PetscMalloc1(7*nvars,&entries);CHKERRQ(ierr);
631f4bdf6c4SBarry Smith   ordering= ex->dofs_order;  /* ordering= 0   nodal ordering
632f4bdf6c4SBarry Smith                                           1   variable ordering */
633f4bdf6c4SBarry Smith   /* stencil entries are orderer by variables: var0_stencil0, var0_stencil1, ..., var0_stencil6, var1_stencil0, var1_stencil1, ...  */
634f4bdf6c4SBarry Smith 
635f4bdf6c4SBarry Smith   if (!ordering) {
636f4bdf6c4SBarry Smith     for (i=0; i<nrow; i++) {
637f4bdf6c4SBarry Smith       grid_rank= irow[i]/nvars;
638f4bdf6c4SBarry Smith       var_type = (irow[i] % nvars);
639f4bdf6c4SBarry Smith 
640f4bdf6c4SBarry Smith       for (j=0; j<ncol; j++) {
641f4bdf6c4SBarry Smith         to_grid_rank= icol[j]/nvars;
642f4bdf6c4SBarry Smith         to_var_type = (icol[j] % nvars);
643f4bdf6c4SBarry Smith 
644f4bdf6c4SBarry Smith         to_var_entry= to_var_entry*7;
645f4bdf6c4SBarry Smith         entries[j]  = to_var_entry;
646f4bdf6c4SBarry Smith 
647f4bdf6c4SBarry Smith         stencil = to_grid_rank-grid_rank;
648f4bdf6c4SBarry Smith         if (!stencil) {
649f4bdf6c4SBarry Smith           entries[j] += 3;
650f4bdf6c4SBarry Smith         } else if (stencil == -1) {
651f4bdf6c4SBarry Smith           entries[j] += 2;
652f4bdf6c4SBarry Smith         } else if (stencil == 1) {
653f4bdf6c4SBarry Smith           entries[j] += 4;
654f4bdf6c4SBarry Smith         } else if (stencil == -ex->gnx) {
655f4bdf6c4SBarry Smith           entries[j] += 1;
656f4bdf6c4SBarry Smith         } else if (stencil == ex->gnx) {
657f4bdf6c4SBarry Smith           entries[j] += 5;
658f4bdf6c4SBarry Smith         } else if (stencil == -ex->gnxgny) {
659f4bdf6c4SBarry Smith           entries[j] += 0;
660f4bdf6c4SBarry Smith         } else if (stencil == ex->gnxgny) {
661f4bdf6c4SBarry Smith           entries[j] += 6;
662f4bdf6c4SBarry 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);
663f4bdf6c4SBarry Smith       }
664f4bdf6c4SBarry Smith 
665f4bdf6c4SBarry Smith       row      = ex->gindices[grid_rank] - ex->rstart;
666f4bdf6c4SBarry Smith       index[0] = ex->xs + (row % ex->nx);
667f4bdf6c4SBarry Smith       index[1] = ex->ys + ((row/ex->nx) % ex->ny);
668f4bdf6c4SBarry Smith       index[2] = ex->zs + (row/(ex->nxny));
669f4bdf6c4SBarry Smith 
670f4bdf6c4SBarry Smith       if (addv == ADD_VALUES) {
6718b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructMatrixAddToValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
672f4bdf6c4SBarry Smith       } else {
6738b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
674f4bdf6c4SBarry Smith       }
675f4bdf6c4SBarry Smith       values += ncol;
676f4bdf6c4SBarry Smith     }
677f4bdf6c4SBarry Smith   } else {
678f4bdf6c4SBarry Smith     for (i=0; i<nrow; i++) {
679f4bdf6c4SBarry Smith       var_type = irow[i]/(ex->gnxgnygnz);
680f4bdf6c4SBarry Smith       grid_rank= irow[i] - var_type*(ex->gnxgnygnz);
681f4bdf6c4SBarry Smith 
682f4bdf6c4SBarry Smith       for (j=0; j<ncol; j++) {
683f4bdf6c4SBarry Smith         to_var_type = icol[j]/(ex->gnxgnygnz);
684f4bdf6c4SBarry Smith         to_grid_rank= icol[j] - to_var_type*(ex->gnxgnygnz);
685f4bdf6c4SBarry Smith 
686f4bdf6c4SBarry Smith         to_var_entry= to_var_entry*7;
687f4bdf6c4SBarry Smith         entries[j]  = to_var_entry;
688f4bdf6c4SBarry Smith 
689f4bdf6c4SBarry Smith         stencil = to_grid_rank-grid_rank;
690f4bdf6c4SBarry Smith         if (!stencil) {
691f4bdf6c4SBarry Smith           entries[j] += 3;
692f4bdf6c4SBarry Smith         } else if (stencil == -1) {
693f4bdf6c4SBarry Smith           entries[j] += 2;
694f4bdf6c4SBarry Smith         } else if (stencil == 1) {
695f4bdf6c4SBarry Smith           entries[j] += 4;
696f4bdf6c4SBarry Smith         } else if (stencil == -ex->gnx) {
697f4bdf6c4SBarry Smith           entries[j] += 1;
698f4bdf6c4SBarry Smith         } else if (stencil == ex->gnx) {
699f4bdf6c4SBarry Smith           entries[j] += 5;
700f4bdf6c4SBarry Smith         } else if (stencil == -ex->gnxgny) {
701f4bdf6c4SBarry Smith           entries[j] += 0;
702f4bdf6c4SBarry Smith         } else if (stencil == ex->gnxgny) {
703f4bdf6c4SBarry Smith           entries[j] += 6;
704f4bdf6c4SBarry 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);
705f4bdf6c4SBarry Smith       }
706f4bdf6c4SBarry Smith 
707f4bdf6c4SBarry Smith       row      = ex->gindices[grid_rank] - ex->rstart;
708f4bdf6c4SBarry Smith       index[0] = ex->xs + (row % ex->nx);
709f4bdf6c4SBarry Smith       index[1] = ex->ys + ((row/ex->nx) % ex->ny);
710f4bdf6c4SBarry Smith       index[2] = ex->zs + (row/(ex->nxny));
711f4bdf6c4SBarry Smith 
712f4bdf6c4SBarry Smith       if (addv == ADD_VALUES) {
7138b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructMatrixAddToValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
714f4bdf6c4SBarry Smith       } else {
7158b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
716f4bdf6c4SBarry Smith       }
717f4bdf6c4SBarry Smith       values += ncol;
718f4bdf6c4SBarry Smith     }
719f4bdf6c4SBarry Smith   }
720f4bdf6c4SBarry Smith   ierr = PetscFree(entries);CHKERRQ(ierr);
721f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
722f4bdf6c4SBarry Smith }
723f4bdf6c4SBarry Smith 
724f4bdf6c4SBarry Smith #undef __FUNCT__
725f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroRowsLocal_HYPRESStruct_3d"
726f4bdf6c4SBarry Smith PetscErrorCode  MatZeroRowsLocal_HYPRESStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscScalar d,Vec x,Vec b)
727f4bdf6c4SBarry Smith {
728f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
729f4bdf6c4SBarry Smith   PetscInt         i,index[3];
730f4bdf6c4SBarry Smith   PetscScalar      **values;
731f4bdf6c4SBarry Smith   Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data;
732f4bdf6c4SBarry Smith 
7334ddd07fcSJed Brown   PetscInt part     = 0;      /* Petsc sstruct interface only allows 1 part */
7344ddd07fcSJed Brown   PetscInt ordering = ex->dofs_order;
7354ddd07fcSJed Brown   PetscInt grid_rank;
7364ddd07fcSJed Brown   PetscInt var_type;
7374ddd07fcSJed Brown   PetscInt nvars= ex->nvars;
738f4bdf6c4SBarry Smith   PetscInt row,*entries;
739f4bdf6c4SBarry Smith 
740f4bdf6c4SBarry Smith   PetscFunctionBegin;
741ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No support");
742785e854fSJed Brown   ierr = PetscMalloc1(7*nvars,&entries);CHKERRQ(ierr);
743f4bdf6c4SBarry Smith 
744785e854fSJed Brown   ierr = PetscMalloc1(nvars,&values);CHKERRQ(ierr);
745785e854fSJed Brown   ierr = PetscMalloc1(7*nvars*nvars,&values[0]);CHKERRQ(ierr);
746f4bdf6c4SBarry Smith   for (i=1; i<nvars; i++) {
747f4bdf6c4SBarry Smith     values[i] = values[i-1] + nvars*7;
748f4bdf6c4SBarry Smith   }
749f4bdf6c4SBarry Smith 
750f4bdf6c4SBarry Smith   for (i=0; i< nvars; i++) {
751f4bdf6c4SBarry Smith     ierr           = PetscMemzero(values[i],nvars*7*sizeof(PetscScalar));CHKERRQ(ierr);
752f4bdf6c4SBarry Smith     *(values[i]+3) = d;
753f4bdf6c4SBarry Smith   }
754f4bdf6c4SBarry Smith 
7558865f1eaSKarl Rupp   for (i= 0; i< nvars*7; i++) entries[i] = i;
756f4bdf6c4SBarry Smith 
757f4bdf6c4SBarry Smith   if (!ordering) {
758f4bdf6c4SBarry Smith     for (i=0; i<nrow; i++) {
759f4bdf6c4SBarry Smith       grid_rank = irow[i] / nvars;
760f4bdf6c4SBarry Smith       var_type  =(irow[i] % nvars);
761f4bdf6c4SBarry Smith 
762f4bdf6c4SBarry Smith       row      = ex->gindices[grid_rank] - ex->rstart;
763f4bdf6c4SBarry Smith       index[0] = ex->xs + (row % ex->nx);
764f4bdf6c4SBarry Smith       index[1] = ex->ys + ((row/ex->nx) % ex->ny);
765f4bdf6c4SBarry Smith       index[2] = ex->zs + (row/(ex->nxny));
7668b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,7*nvars,(HYPRE_Int *)entries,values[var_type]));
767f4bdf6c4SBarry Smith     }
768f4bdf6c4SBarry Smith   } else {
769f4bdf6c4SBarry Smith     for (i=0; i<nrow; i++) {
770f4bdf6c4SBarry Smith       var_type = irow[i]/(ex->gnxgnygnz);
771f4bdf6c4SBarry Smith       grid_rank= irow[i] - var_type*(ex->gnxgnygnz);
772f4bdf6c4SBarry Smith 
773f4bdf6c4SBarry Smith       row      = ex->gindices[grid_rank] - ex->rstart;
774f4bdf6c4SBarry Smith       index[0] = ex->xs + (row % ex->nx);
775f4bdf6c4SBarry Smith       index[1] = ex->ys + ((row/ex->nx) % ex->ny);
776f4bdf6c4SBarry Smith       index[2] = ex->zs + (row/(ex->nxny));
7778b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,7*nvars,(HYPRE_Int *)entries,values[var_type]));
778f4bdf6c4SBarry Smith     }
779f4bdf6c4SBarry Smith   }
780fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructMatrixAssemble,(ex->ss_mat));
781f4bdf6c4SBarry Smith   ierr = PetscFree(values[0]);CHKERRQ(ierr);
782f4bdf6c4SBarry Smith   ierr = PetscFree(values);CHKERRQ(ierr);
783f4bdf6c4SBarry Smith   ierr = PetscFree(entries);CHKERRQ(ierr);
784f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
785f4bdf6c4SBarry Smith }
786f4bdf6c4SBarry Smith 
787f4bdf6c4SBarry Smith #undef __FUNCT__
788f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPRESStruct_3d"
789f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPRESStruct_3d(Mat mat)
790f4bdf6c4SBarry Smith {
791f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
792f4bdf6c4SBarry Smith   Mat_HYPRESStruct *ex  = (Mat_HYPRESStruct*) mat->data;
7934ddd07fcSJed Brown   PetscInt         nvars= ex->nvars;
7944ddd07fcSJed Brown   PetscInt         size;
7954ddd07fcSJed Brown   PetscInt         part= 0;   /* only one part */
796f4bdf6c4SBarry Smith 
797f4bdf6c4SBarry Smith   PetscFunctionBegin;
798f4bdf6c4SBarry 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);
799f4bdf6c4SBarry Smith   {
8004ddd07fcSJed Brown     PetscInt    i,*entries,iupper[3],ilower[3];
801f4bdf6c4SBarry Smith     PetscScalar *values;
8024ddd07fcSJed Brown 
803f4bdf6c4SBarry Smith 
804f4bdf6c4SBarry Smith     for (i= 0; i< 3; i++) {
805f4bdf6c4SBarry Smith       ilower[i]= ex->hbox.imin[i];
806f4bdf6c4SBarry Smith       iupper[i]= ex->hbox.imax[i];
807f4bdf6c4SBarry Smith     }
808f4bdf6c4SBarry Smith 
809dcca6d9dSJed Brown     ierr = PetscMalloc2(nvars*7,&entries,nvars*7*size,&values);CHKERRQ(ierr);
8108865f1eaSKarl Rupp     for (i= 0; i< nvars*7; i++) entries[i]= i;
811f4bdf6c4SBarry Smith     ierr = PetscMemzero(values,nvars*7*size*sizeof(PetscScalar));CHKERRQ(ierr);
812f4bdf6c4SBarry Smith 
813f4bdf6c4SBarry Smith     for (i= 0; i< nvars; i++) {
8148b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructMatrixSetBoxValues,(ex->ss_mat,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,nvars*7,(HYPRE_Int *)entries,values));
815f4bdf6c4SBarry Smith     }
816f4bdf6c4SBarry Smith     ierr = PetscFree2(entries,values);CHKERRQ(ierr);
817f4bdf6c4SBarry Smith   }
818fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructMatrixAssemble,(ex->ss_mat));
819f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
820f4bdf6c4SBarry Smith }
821f4bdf6c4SBarry Smith 
822f4bdf6c4SBarry Smith 
823f4bdf6c4SBarry Smith #undef __FUNCT__
824c688c046SMatthew G Knepley #define __FUNCT__ "MatSetupDM_HYPRESStruct"
8251e6b0712SBarry Smith static PetscErrorCode  MatSetupDM_HYPRESStruct(Mat mat,DM da)
826f4bdf6c4SBarry Smith {
827f4bdf6c4SBarry Smith   PetscErrorCode         ierr;
828f4bdf6c4SBarry Smith   Mat_HYPRESStruct       *ex = (Mat_HYPRESStruct*) mat->data;
829f4bdf6c4SBarry Smith   PetscInt               dim,dof,sw[3],nx,ny,nz;
8304ddd07fcSJed Brown   PetscInt               ilower[3],iupper[3],ssize,i;
831bff4a2f0SMatthew G. Knepley   DMBoundaryType         px,py,pz;
832f4bdf6c4SBarry Smith   DMDAStencilType        st;
8334ddd07fcSJed Brown   PetscInt               nparts= 1;  /* assuming only one part */
8344ddd07fcSJed Brown   PetscInt               part  = 0;
835565245c5SBarry Smith   ISLocalToGlobalMapping ltog;
836f4bdf6c4SBarry Smith   PetscFunctionBegin;
837f4bdf6c4SBarry Smith   ex->da = da;
838f4bdf6c4SBarry Smith   ierr   = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);
839f4bdf6c4SBarry Smith 
8407ab44359SJed Brown   ierr       = DMDAGetInfo(ex->da,&dim,0,0,0,0,0,0,&dof,&sw[0],&px,&py,&pz,&st);CHKERRQ(ierr);
841f4bdf6c4SBarry Smith   ierr       = DMDAGetCorners(ex->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
842f4bdf6c4SBarry Smith   iupper[0] += ilower[0] - 1;
843f4bdf6c4SBarry Smith   iupper[1] += ilower[1] - 1;
844f4bdf6c4SBarry Smith   iupper[2] += ilower[2] - 1;
845f4bdf6c4SBarry Smith   /* the hypre_Box is used to zero out the matrix entries in MatZeroValues() */
846f4bdf6c4SBarry Smith   ex->hbox.imin[0] = ilower[0];
847f4bdf6c4SBarry Smith   ex->hbox.imin[1] = ilower[1];
848f4bdf6c4SBarry Smith   ex->hbox.imin[2] = ilower[2];
849f4bdf6c4SBarry Smith   ex->hbox.imax[0] = iupper[0];
850f4bdf6c4SBarry Smith   ex->hbox.imax[1] = iupper[1];
851f4bdf6c4SBarry Smith   ex->hbox.imax[2] = iupper[2];
852f4bdf6c4SBarry Smith 
853f4bdf6c4SBarry Smith   ex->dofs_order = 0;
854f4bdf6c4SBarry Smith 
855f4bdf6c4SBarry Smith   /* assuming that the same number of dofs on each gridpoint. Also assume all cell-centred based */
856f4bdf6c4SBarry Smith   ex->nvars= dof;
857f4bdf6c4SBarry Smith 
858f4bdf6c4SBarry Smith   /* create the hypre grid object and set its information */
859ce94432eSBarry Smith   if (px || py || pz) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add periodic support by calling HYPRE_SStructGridSetPeriodic()");
860fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGridCreate,(ex->hcomm,dim,nparts,&ex->ss_grid));
861f4bdf6c4SBarry Smith 
862fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGridSetExtents,(ex->ss_grid,part,ex->hbox.imin,ex->hbox.imax));
863f4bdf6c4SBarry Smith 
864f4bdf6c4SBarry Smith   {
865f4bdf6c4SBarry Smith     HYPRE_SStructVariable *vartypes;
866785e854fSJed Brown     ierr = PetscMalloc1(ex->nvars,&vartypes);CHKERRQ(ierr);
8678865f1eaSKarl Rupp     for (i= 0; i< ex->nvars; i++) vartypes[i]= HYPRE_SSTRUCT_VARIABLE_CELL;
868fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructGridSetVariables,(ex->ss_grid, part, ex->nvars,vartypes));
869f4bdf6c4SBarry Smith     ierr = PetscFree(vartypes);CHKERRQ(ierr);
870f4bdf6c4SBarry Smith   }
871fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGridAssemble,(ex->ss_grid));
872f4bdf6c4SBarry Smith 
873f4bdf6c4SBarry Smith   sw[1] = sw[0];
874f4bdf6c4SBarry Smith   sw[2] = sw[1];
875fd3f9acdSBarry Smith   /* PetscStackCallStandard(HYPRE_SStructGridSetNumGhost,(ex->ss_grid,sw)); */
876f4bdf6c4SBarry Smith 
877f4bdf6c4SBarry Smith   /* create the hypre stencil object and set its information */
878ce94432eSBarry Smith   if (sw[0] > 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add support for wider stencils");
879ce94432eSBarry Smith   if (st == DMDA_STENCIL_BOX) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add support for box stencils");
880f4bdf6c4SBarry Smith 
881f4bdf6c4SBarry Smith   if (dim == 1) {
8824ddd07fcSJed Brown     PetscInt offsets[3][1] = {{-1},{0},{1}};
8834ddd07fcSJed Brown     PetscInt j, cnt;
884f4bdf6c4SBarry Smith 
885f4bdf6c4SBarry Smith     ssize = 3*(ex->nvars);
886fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil));
887f4bdf6c4SBarry Smith     cnt= 0;
888f4bdf6c4SBarry Smith     for (i = 0; i < (ex->nvars); i++) {
889f4bdf6c4SBarry Smith       for (j = 0; j < 3; j++) {
8908b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, (HYPRE_Int *)offsets[j], i));
891f4bdf6c4SBarry Smith         cnt++;
892f4bdf6c4SBarry Smith       }
893f4bdf6c4SBarry Smith     }
894f4bdf6c4SBarry Smith   } else if (dim == 2) {
8954ddd07fcSJed Brown     PetscInt offsets[5][2] = {{0,-1},{-1,0},{0,0},{1,0},{0,1}};
8964ddd07fcSJed Brown     PetscInt j, cnt;
897f4bdf6c4SBarry Smith 
898f4bdf6c4SBarry Smith     ssize = 5*(ex->nvars);
899fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil));
900f4bdf6c4SBarry Smith     cnt= 0;
901f4bdf6c4SBarry Smith     for (i= 0; i< (ex->nvars); i++) {
902f4bdf6c4SBarry Smith       for (j= 0; j< 5; j++) {
9038b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, (HYPRE_Int *)offsets[j], i));
904f4bdf6c4SBarry Smith         cnt++;
905f4bdf6c4SBarry Smith       }
906f4bdf6c4SBarry Smith     }
907f4bdf6c4SBarry Smith   } else if (dim == 3) {
9084ddd07fcSJed Brown     PetscInt offsets[7][3] = {{0,0,-1},{0,-1,0},{-1,0,0},{0,0,0},{1,0,0},{0,1,0},{0,0,1}};
9094ddd07fcSJed Brown     PetscInt j, cnt;
910f4bdf6c4SBarry Smith 
911f4bdf6c4SBarry Smith     ssize = 7*(ex->nvars);
912fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil));
913f4bdf6c4SBarry Smith     cnt= 0;
914f4bdf6c4SBarry Smith     for (i= 0; i< (ex->nvars); i++) {
915f4bdf6c4SBarry Smith       for (j= 0; j< 7; j++) {
9168b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, (HYPRE_Int *)offsets[j], i));
917f4bdf6c4SBarry Smith         cnt++;
918f4bdf6c4SBarry Smith       }
919f4bdf6c4SBarry Smith     }
920f4bdf6c4SBarry Smith   }
921f4bdf6c4SBarry Smith 
922f4bdf6c4SBarry Smith   /* create the HYPRE graph */
923fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGraphCreate,(ex->hcomm, ex->ss_grid, &(ex->ss_graph)));
924f4bdf6c4SBarry Smith 
925f4bdf6c4SBarry Smith   /* set the stencil graph. Note that each variable has the same graph. This means that each
926f4bdf6c4SBarry Smith      variable couples to all the other variable and with the same stencil pattern. */
927f4bdf6c4SBarry Smith   for (i= 0; i< (ex->nvars); i++) {
928fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructGraphSetStencil,(ex->ss_graph,part,i,ex->ss_stencil));
929f4bdf6c4SBarry Smith   }
930fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGraphAssemble,(ex->ss_graph));
931f4bdf6c4SBarry Smith 
932f4bdf6c4SBarry Smith   /* create the HYPRE sstruct vectors for rhs and solution */
933fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorCreate,(ex->hcomm,ex->ss_grid,&ex->ss_b));
934fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorCreate,(ex->hcomm,ex->ss_grid,&ex->ss_x));
935fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorInitialize,(ex->ss_b));
936fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorInitialize,(ex->ss_x));
937fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorAssemble,(ex->ss_b));
938fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorAssemble,(ex->ss_x));
939f4bdf6c4SBarry Smith 
940f4bdf6c4SBarry Smith   /* create the hypre matrix object and set its information */
941fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructMatrixCreate,(ex->hcomm,ex->ss_graph,&ex->ss_mat));
942fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGridDestroy,(ex->ss_grid));
943fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructStencilDestroy,(ex->ss_stencil));
944f4bdf6c4SBarry Smith   if (ex->needsinitialization) {
945fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructMatrixInitialize,(ex->ss_mat));
946f4bdf6c4SBarry Smith     ex->needsinitialization = PETSC_FALSE;
947f4bdf6c4SBarry Smith   }
948f4bdf6c4SBarry Smith 
949f4bdf6c4SBarry Smith   /* set the global and local sizes of the matrix */
950f4bdf6c4SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
951f4bdf6c4SBarry Smith   ierr = MatSetSizes(mat,dof*nx*ny*nz,dof*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
952f4bdf6c4SBarry Smith   ierr = PetscLayoutSetBlockSize(mat->rmap,1);CHKERRQ(ierr);
953f4bdf6c4SBarry Smith   ierr = PetscLayoutSetBlockSize(mat->cmap,1);CHKERRQ(ierr);
954f4bdf6c4SBarry Smith   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
955f4bdf6c4SBarry Smith   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
956f4bdf6c4SBarry Smith 
957f4bdf6c4SBarry Smith   if (dim == 3) {
958f4bdf6c4SBarry Smith     mat->ops->setvalueslocal = MatSetValuesLocal_HYPRESStruct_3d;
959f4bdf6c4SBarry Smith     mat->ops->zerorowslocal  = MatZeroRowsLocal_HYPRESStruct_3d;
960f4bdf6c4SBarry Smith     mat->ops->zeroentries    = MatZeroEntries_HYPRESStruct_3d;
9618865f1eaSKarl Rupp 
962f4bdf6c4SBarry Smith     ierr = MatZeroEntries_HYPRESStruct_3d(mat);CHKERRQ(ierr);
963ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only support for 3d DMDA currently");
964f4bdf6c4SBarry Smith 
965f4bdf6c4SBarry Smith   /* get values that will be used repeatedly in MatSetValuesLocal() and MatZeroRowsLocal() repeatedly */
9660298fd71SBarry Smith   ierr = MatGetOwnershipRange(mat,&ex->rstart,NULL);CHKERRQ(ierr);
967565245c5SBarry Smith   ierr = DMGetLocalToGlobalMapping(ex->da,&ltog);CHKERRQ(ierr);
968565245c5SBarry Smith   ierr = ISLocalToGlobalMappingGetIndices(ltog, (const PetscInt **) &ex->gindices);CHKERRQ(ierr);
969f4bdf6c4SBarry Smith   ierr = DMDAGetGhostCorners(ex->da,0,0,0,&ex->gnx,&ex->gnxgny,&ex->gnxgnygnz);CHKERRQ(ierr);
9708865f1eaSKarl Rupp 
971f4bdf6c4SBarry Smith   ex->gnxgny    *= ex->gnx;
972f4bdf6c4SBarry Smith   ex->gnxgnygnz *= ex->gnxgny;
9738865f1eaSKarl Rupp 
974f4bdf6c4SBarry Smith   ierr = DMDAGetCorners(ex->da,&ex->xs,&ex->ys,&ex->zs,&ex->nx,&ex->ny,&ex->nz);CHKERRQ(ierr);
9758865f1eaSKarl Rupp 
976f4bdf6c4SBarry Smith   ex->nxny   = ex->nx*ex->ny;
977f4bdf6c4SBarry Smith   ex->nxnynz = ex->nz*ex->nxny;
978f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
979f4bdf6c4SBarry Smith }
980f4bdf6c4SBarry Smith 
981f4bdf6c4SBarry Smith #undef __FUNCT__
982f4bdf6c4SBarry Smith #define __FUNCT__ "MatMult_HYPRESStruct"
983f4bdf6c4SBarry Smith PetscErrorCode MatMult_HYPRESStruct(Mat A,Vec x,Vec y)
984f4bdf6c4SBarry Smith {
985f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
986f4bdf6c4SBarry Smith   PetscScalar      *xx,*yy;
9874ddd07fcSJed Brown   PetscInt         ilower[3],iupper[3];
988f4bdf6c4SBarry Smith   Mat_HYPRESStruct *mx     = (Mat_HYPRESStruct*)(A->data);
9894ddd07fcSJed Brown   PetscInt         ordering= mx->dofs_order;
9904ddd07fcSJed Brown   PetscInt         nvars   = mx->nvars;
9914ddd07fcSJed Brown   PetscInt         part    = 0;
9924ddd07fcSJed Brown   PetscInt         size;
9934ddd07fcSJed Brown   PetscInt         i;
994f4bdf6c4SBarry Smith 
995f4bdf6c4SBarry Smith   PetscFunctionBegin;
996f4bdf6c4SBarry Smith   ierr       = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
997f4bdf6c4SBarry Smith   iupper[0] += ilower[0] - 1;
998f4bdf6c4SBarry Smith   iupper[1] += ilower[1] - 1;
999f4bdf6c4SBarry Smith   iupper[2] += ilower[2] - 1;
1000f4bdf6c4SBarry Smith 
1001f4bdf6c4SBarry Smith   size= 1;
10028865f1eaSKarl Rupp   for (i = 0; i < 3; i++) size *= (iupper[i]-ilower[i]+1);
1003f4bdf6c4SBarry Smith 
1004f4bdf6c4SBarry Smith   /* copy x values over to hypre for variable ordering */
1005f4bdf6c4SBarry Smith   if (ordering) {
1006fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1007f4bdf6c4SBarry Smith     ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1008f4bdf6c4SBarry Smith     for (i= 0; i< nvars; i++) {
10098b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,xx+(size*i)));
1010f4bdf6c4SBarry Smith     }
1011f4bdf6c4SBarry Smith     ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1012fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
1013fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
1014f4bdf6c4SBarry Smith 
1015f4bdf6c4SBarry Smith     /* copy solution values back to PETSc */
1016f4bdf6c4SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1017f4bdf6c4SBarry Smith     for (i= 0; i< nvars; i++) {
10188b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,yy+(size*i)));
1019f4bdf6c4SBarry Smith     }
1020f4bdf6c4SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1021f4bdf6c4SBarry Smith   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
1022f4bdf6c4SBarry Smith     PetscScalar *z;
10234ddd07fcSJed Brown     PetscInt    j, k;
1024f4bdf6c4SBarry Smith 
1025785e854fSJed Brown     ierr = PetscMalloc1(nvars*size,&z);CHKERRQ(ierr);
1026fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
1027f4bdf6c4SBarry Smith     ierr = VecGetArray(x,&xx);CHKERRQ(ierr);
1028f4bdf6c4SBarry Smith 
1029f4bdf6c4SBarry Smith     /* transform nodal to hypre's variable ordering for sys_pfmg */
1030f4bdf6c4SBarry Smith     for (i= 0; i< size; i++) {
1031f4bdf6c4SBarry Smith       k= i*nvars;
10328865f1eaSKarl Rupp       for (j= 0; j< nvars; j++) z[j*size+i]= xx[k+j];
1033f4bdf6c4SBarry Smith     }
1034f4bdf6c4SBarry Smith     for (i= 0; i< nvars; i++) {
10358b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
1036f4bdf6c4SBarry Smith     }
1037f4bdf6c4SBarry Smith     ierr = VecRestoreArray(x,&xx);CHKERRQ(ierr);
1038fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
1039fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
1040f4bdf6c4SBarry Smith 
1041f4bdf6c4SBarry Smith     /* copy solution values back to PETSc */
1042f4bdf6c4SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
1043f4bdf6c4SBarry Smith     for (i= 0; i< nvars; i++) {
10448b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
1045f4bdf6c4SBarry Smith     }
1046f4bdf6c4SBarry Smith     /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
1047f4bdf6c4SBarry Smith     for (i= 0; i< size; i++) {
1048f4bdf6c4SBarry Smith       k= i*nvars;
10498865f1eaSKarl Rupp       for (j= 0; j< nvars; j++) yy[k+j]= z[j*size+i];
1050f4bdf6c4SBarry Smith     }
1051f4bdf6c4SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
1052f4bdf6c4SBarry Smith     ierr = PetscFree(z);CHKERRQ(ierr);
1053f4bdf6c4SBarry Smith   }
1054f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
1055f4bdf6c4SBarry Smith }
1056f4bdf6c4SBarry Smith 
1057f4bdf6c4SBarry Smith #undef __FUNCT__
1058f4bdf6c4SBarry Smith #define __FUNCT__ "MatAssemblyEnd_HYPRESStruct"
1059f4bdf6c4SBarry Smith PetscErrorCode MatAssemblyEnd_HYPRESStruct(Mat mat,MatAssemblyType mode)
1060f4bdf6c4SBarry Smith {
1061f4bdf6c4SBarry Smith   Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data;
1062f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
1063f4bdf6c4SBarry Smith 
1064f4bdf6c4SBarry Smith   PetscFunctionBegin;
1065fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructMatrixAssemble,(ex->ss_mat));
1066f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
1067f4bdf6c4SBarry Smith }
1068f4bdf6c4SBarry Smith 
1069f4bdf6c4SBarry Smith #undef __FUNCT__
1070f4bdf6c4SBarry Smith #define __FUNCT__ "MatZeroEntries_HYPRESStruct"
1071f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPRESStruct(Mat mat)
1072f4bdf6c4SBarry Smith {
1073f4bdf6c4SBarry Smith   PetscFunctionBegin;
1074f4bdf6c4SBarry Smith   /* before the DMDA is set to the matrix the zero doesn't need to do anything */
1075f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
1076f4bdf6c4SBarry Smith }
1077f4bdf6c4SBarry Smith 
1078f4bdf6c4SBarry Smith 
1079f4bdf6c4SBarry Smith #undef __FUNCT__
1080f4bdf6c4SBarry Smith #define __FUNCT__ "MatDestroy_HYPRESStruct"
1081f4bdf6c4SBarry Smith PetscErrorCode MatDestroy_HYPRESStruct(Mat mat)
1082f4bdf6c4SBarry Smith {
1083f4bdf6c4SBarry Smith   Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data;
1084f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
1085f4bdf6c4SBarry Smith 
1086f4bdf6c4SBarry Smith   PetscFunctionBegin;
1087fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGraphDestroy,(ex->ss_graph));
1088fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructMatrixDestroy,(ex->ss_mat));
1089fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorDestroy,(ex->ss_x));
1090fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorDestroy,(ex->ss_b));
1091f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
1092f4bdf6c4SBarry Smith }
1093f4bdf6c4SBarry Smith 
1094f4bdf6c4SBarry Smith #undef __FUNCT__
1095f4bdf6c4SBarry Smith #define __FUNCT__ "MatCreate_HYPRESStruct"
10968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_HYPRESStruct(Mat B)
1097f4bdf6c4SBarry Smith {
1098f4bdf6c4SBarry Smith   Mat_HYPRESStruct *ex;
1099f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
1100f4bdf6c4SBarry Smith 
1101f4bdf6c4SBarry Smith   PetscFunctionBegin;
1102b00a9115SJed Brown   ierr         = PetscNewLog(B,&ex);CHKERRQ(ierr);
1103f4bdf6c4SBarry Smith   B->data      = (void*)ex;
1104f4bdf6c4SBarry Smith   B->rmap->bs  = 1;
1105f4bdf6c4SBarry Smith   B->assembled = PETSC_FALSE;
1106f4bdf6c4SBarry Smith 
1107f4bdf6c4SBarry Smith   B->insertmode = NOT_SET_VALUES;
1108f4bdf6c4SBarry Smith 
1109f4bdf6c4SBarry Smith   B->ops->assemblyend = MatAssemblyEnd_HYPRESStruct;
1110f4bdf6c4SBarry Smith   B->ops->mult        = MatMult_HYPRESStruct;
1111f4bdf6c4SBarry Smith   B->ops->zeroentries = MatZeroEntries_HYPRESStruct;
1112f4bdf6c4SBarry Smith   B->ops->destroy     = MatDestroy_HYPRESStruct;
1113f4bdf6c4SBarry Smith 
1114f4bdf6c4SBarry Smith   ex->needsinitialization = PETSC_TRUE;
1115f4bdf6c4SBarry Smith 
1116ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&(ex->hcomm));CHKERRQ(ierr);
1117bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSetupDM_C",MatSetupDM_HYPRESStruct);CHKERRQ(ierr);
1118f4bdf6c4SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRESSTRUCT);CHKERRQ(ierr);
1119f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
1120f4bdf6c4SBarry Smith }
1121f4bdf6c4SBarry Smith 
1122f4bdf6c4SBarry Smith 
1123f4bdf6c4SBarry Smith 
1124