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