xref: /petsc/src/dm/impls/da/hypre/mhyp.c (revision b026d285daece26300d337cd10f27db66bb146fb)
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 /* -----------------------------------------------------------------------------------------------------------------*/
11f4bdf6c4SBarry Smith 
12f4bdf6c4SBarry Smith /*MC
13f4bdf6c4SBarry Smith    MATHYPRESTRUCT - MATHYPRESTRUCT = "hyprestruct" - A matrix type to be used for parallel sparse matrices
14f4bdf6c4SBarry Smith           based on the hypre HYPRE_StructMatrix.
15f4bdf6c4SBarry Smith 
16f4bdf6c4SBarry Smith    Level: intermediate
17f4bdf6c4SBarry Smith 
18f4bdf6c4SBarry Smith    Notes: Unlike the more general support for blocks in hypre this allows only one block per process and requires the block
19f4bdf6c4SBarry Smith           be defined by a DMDA.
20f4bdf6c4SBarry Smith 
21c688c046SMatthew G Knepley           The matrix needs a DMDA associated with it by either a call to MatSetupDM() or if the matrix is obtained from DMCreateMatrix()
22f4bdf6c4SBarry Smith 
23c688c046SMatthew G Knepley .seealso: MatCreate(), PCPFMG, MatSetupDM(), DMCreateMatrix()
24f4bdf6c4SBarry Smith M*/
25f4bdf6c4SBarry Smith 
26f4bdf6c4SBarry Smith 
27f4bdf6c4SBarry Smith PetscErrorCode  MatSetValuesLocal_HYPREStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
28f4bdf6c4SBarry Smith {
29f4bdf6c4SBarry Smith   PetscErrorCode    ierr;
30f4bdf6c4SBarry Smith   PetscInt          i,j,stencil,index[3],row,entries[7];
31f4bdf6c4SBarry Smith   const PetscScalar *values = y;
32f4bdf6c4SBarry Smith   Mat_HYPREStruct   *ex     = (Mat_HYPREStruct*) mat->data;
33f4bdf6c4SBarry Smith 
34f4bdf6c4SBarry Smith   PetscFunctionBegin;
35f6dae198SJed Brown   if (PetscUnlikely(ncol >= 7)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"ncol %D >= 7 too large",ncol);
36f4bdf6c4SBarry Smith   for (i=0; i<nrow; i++) {
37f4bdf6c4SBarry Smith     for (j=0; j<ncol; j++) {
38f4bdf6c4SBarry Smith       stencil = icol[j] - irow[i];
39f4bdf6c4SBarry Smith       if (!stencil) {
40f4bdf6c4SBarry Smith         entries[j] = 3;
41f4bdf6c4SBarry Smith       } else if (stencil == -1) {
42f4bdf6c4SBarry Smith         entries[j] = 2;
43f4bdf6c4SBarry Smith       } else if (stencil == 1) {
44f4bdf6c4SBarry Smith         entries[j] = 4;
45f4bdf6c4SBarry Smith       } else if (stencil == -ex->gnx) {
46f4bdf6c4SBarry Smith         entries[j] = 1;
47f4bdf6c4SBarry Smith       } else if (stencil == ex->gnx) {
48f4bdf6c4SBarry Smith         entries[j] = 5;
49f4bdf6c4SBarry Smith       } else if (stencil == -ex->gnxgny) {
50f4bdf6c4SBarry Smith         entries[j] = 0;
51f4bdf6c4SBarry Smith       } else if (stencil == ex->gnxgny) {
52f4bdf6c4SBarry Smith         entries[j] = 6;
53f4bdf6c4SBarry 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);
54f4bdf6c4SBarry Smith     }
55f4bdf6c4SBarry Smith     row      = ex->gindices[irow[i]] - ex->rstart;
56f4bdf6c4SBarry Smith     index[0] = ex->xs + (row % ex->nx);
57f4bdf6c4SBarry Smith     index[1] = ex->ys + ((row/ex->nx) % ex->ny);
58f4bdf6c4SBarry Smith     index[2] = ex->zs + (row/(ex->nxny));
59f4bdf6c4SBarry Smith     if (addv == ADD_VALUES) {
608b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_StructMatrixAddToValues,(ex->hmat,(HYPRE_Int *)index,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
61f4bdf6c4SBarry Smith     } else {
628b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_StructMatrixSetValues,(ex->hmat,(HYPRE_Int *)index,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
63f4bdf6c4SBarry Smith     }
64f4bdf6c4SBarry Smith     values += ncol;
65f4bdf6c4SBarry Smith   }
66f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
67f4bdf6c4SBarry Smith }
68f4bdf6c4SBarry Smith 
69f4bdf6c4SBarry Smith PetscErrorCode  MatZeroRowsLocal_HYPREStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscScalar d,Vec x,Vec b)
70f4bdf6c4SBarry Smith {
71f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
72f4bdf6c4SBarry Smith   PetscInt        i,index[3],row,entries[7] = {0,1,2,3,4,5,6};
73f4bdf6c4SBarry Smith   PetscScalar     values[7];
74f4bdf6c4SBarry Smith   Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data;
75f4bdf6c4SBarry Smith 
76f4bdf6c4SBarry Smith   PetscFunctionBegin;
77ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No support");
78f4bdf6c4SBarry Smith   ierr      = PetscMemzero(values,7*sizeof(PetscScalar));CHKERRQ(ierr);
79f4bdf6c4SBarry Smith   values[3] = d;
80f4bdf6c4SBarry Smith   for (i=0; i<nrow; i++) {
81f4bdf6c4SBarry Smith     row      = ex->gindices[irow[i]] - ex->rstart;
82f4bdf6c4SBarry Smith     index[0] = ex->xs + (row % ex->nx);
83f4bdf6c4SBarry Smith     index[1] = ex->ys + ((row/ex->nx) % ex->ny);
84f4bdf6c4SBarry Smith     index[2] = ex->zs + (row/(ex->nxny));
858b1f7689SBarry Smith     PetscStackCallStandard(HYPRE_StructMatrixSetValues,(ex->hmat,(HYPRE_Int *)index,7,(HYPRE_Int *)entries,values));
86f4bdf6c4SBarry Smith   }
87fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixAssemble,(ex->hmat));
88f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
89f4bdf6c4SBarry Smith }
90f4bdf6c4SBarry Smith 
91f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPREStruct_3d(Mat mat)
92f4bdf6c4SBarry Smith {
93f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
94f4bdf6c4SBarry Smith   PetscInt        indices[7] = {0,1,2,3,4,5,6};
95f4bdf6c4SBarry Smith   Mat_HYPREStruct *ex        = (Mat_HYPREStruct*) mat->data;
96f4bdf6c4SBarry Smith 
97f4bdf6c4SBarry Smith   PetscFunctionBegin;
98f4bdf6c4SBarry Smith   /* hypre has no public interface to do this */
998b1f7689SBarry Smith   PetscStackCallStandard(hypre_StructMatrixClearBoxValues,(ex->hmat,&ex->hbox,7,(HYPRE_Int *)indices,0,1));
100fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixAssemble,(ex->hmat));
101f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
102f4bdf6c4SBarry Smith }
103f4bdf6c4SBarry Smith 
1041e6b0712SBarry Smith static PetscErrorCode  MatSetupDM_HYPREStruct(Mat mat,DM da)
105f4bdf6c4SBarry Smith {
106f4bdf6c4SBarry Smith   PetscErrorCode         ierr;
107f4bdf6c4SBarry Smith   Mat_HYPREStruct        *ex = (Mat_HYPREStruct*) mat->data;
1084ddd07fcSJed Brown   PetscInt               dim,dof,sw[3],nx,ny,nz,ilower[3],iupper[3],ssize,i;
109bff4a2f0SMatthew G. Knepley   DMBoundaryType         px,py,pz;
110f4bdf6c4SBarry Smith   DMDAStencilType        st;
111565245c5SBarry Smith   ISLocalToGlobalMapping ltog;
112f4bdf6c4SBarry Smith 
113f4bdf6c4SBarry Smith   PetscFunctionBegin;
114f4bdf6c4SBarry Smith   ex->da = da;
115f4bdf6c4SBarry Smith   ierr   = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);
116f4bdf6c4SBarry Smith 
1177ab44359SJed Brown   ierr       = DMDAGetInfo(ex->da,&dim,0,0,0,0,0,0,&dof,&sw[0],&px,&py,&pz,&st);CHKERRQ(ierr);
118f4bdf6c4SBarry Smith   ierr       = DMDAGetCorners(ex->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
119f4bdf6c4SBarry Smith   iupper[0] += ilower[0] - 1;
120f4bdf6c4SBarry Smith   iupper[1] += ilower[1] - 1;
121f4bdf6c4SBarry Smith   iupper[2] += ilower[2] - 1;
122f4bdf6c4SBarry Smith 
123f4bdf6c4SBarry Smith   /* the hypre_Box is used to zero out the matrix entries in MatZeroValues() */
124f4bdf6c4SBarry Smith   ex->hbox.imin[0] = ilower[0];
125f4bdf6c4SBarry Smith   ex->hbox.imin[1] = ilower[1];
126f4bdf6c4SBarry Smith   ex->hbox.imin[2] = ilower[2];
127f4bdf6c4SBarry Smith   ex->hbox.imax[0] = iupper[0];
128f4bdf6c4SBarry Smith   ex->hbox.imax[1] = iupper[1];
129f4bdf6c4SBarry Smith   ex->hbox.imax[2] = iupper[2];
130f4bdf6c4SBarry Smith 
131f4bdf6c4SBarry Smith   /* create the hypre grid object and set its information */
132ce94432eSBarry Smith   if (dof > 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Currently only support for scalar problems");
133ce94432eSBarry Smith   if (px || py || pz) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add periodic support by calling HYPRE_StructGridSetPeriodic()");
134fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructGridCreate,(ex->hcomm,dim,&ex->hgrid));
1358b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_StructGridSetExtents,(ex->hgrid,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper));
136fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructGridAssemble,(ex->hgrid));
137f4bdf6c4SBarry Smith 
138f4bdf6c4SBarry Smith   sw[1] = sw[0];
139f4bdf6c4SBarry Smith   sw[2] = sw[1];
1408b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_StructGridSetNumGhost,(ex->hgrid,(HYPRE_Int *)sw));
141f4bdf6c4SBarry Smith 
142f4bdf6c4SBarry Smith   /* create the hypre stencil object and set its information */
143ce94432eSBarry Smith   if (sw[0] > 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add support for wider stencils");
144ce94432eSBarry Smith   if (st == DMDA_STENCIL_BOX) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add support for box stencils");
145f4bdf6c4SBarry Smith   if (dim == 1) {
1464ddd07fcSJed Brown     PetscInt offsets[3][1] = {{-1},{0},{1}};
147f4bdf6c4SBarry Smith     ssize = 3;
148fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil));
149f4bdf6c4SBarry Smith     for (i=0; i<ssize; i++) {
1508b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_StructStencilSetElement,(ex->hstencil,i,(HYPRE_Int *)offsets[i]));
151f4bdf6c4SBarry Smith     }
152f4bdf6c4SBarry Smith   } else if (dim == 2) {
1534ddd07fcSJed Brown     PetscInt offsets[5][2] = {{0,-1},{-1,0},{0,0},{1,0},{0,1}};
154f4bdf6c4SBarry Smith     ssize = 5;
155fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil));
156f4bdf6c4SBarry Smith     for (i=0; i<ssize; i++) {
1578b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_StructStencilSetElement,(ex->hstencil,i,(HYPRE_Int *)offsets[i]));
158f4bdf6c4SBarry Smith     }
159f4bdf6c4SBarry Smith   } else if (dim == 3) {
1604ddd07fcSJed 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}};
161f4bdf6c4SBarry Smith     ssize = 7;
162fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructStencilCreate,(dim,ssize,&ex->hstencil));
163f4bdf6c4SBarry Smith     for (i=0; i<ssize; i++) {
1648b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_StructStencilSetElement,(ex->hstencil,i,(HYPRE_Int *)offsets[i]));
165f4bdf6c4SBarry Smith     }
166f4bdf6c4SBarry Smith   }
167f4bdf6c4SBarry Smith 
168f4bdf6c4SBarry Smith   /* create the HYPRE vector for rhs and solution */
169fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorCreate,(ex->hcomm,ex->hgrid,&ex->hb));
170fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorCreate,(ex->hcomm,ex->hgrid,&ex->hx));
171fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorInitialize,(ex->hb));
172fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorInitialize,(ex->hx));
173fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorAssemble,(ex->hb));
174fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorAssemble,(ex->hx));
175f4bdf6c4SBarry Smith 
176f4bdf6c4SBarry Smith   /* create the hypre matrix object and set its information */
177fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixCreate,(ex->hcomm,ex->hgrid,ex->hstencil,&ex->hmat));
178fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructGridDestroy,(ex->hgrid));
179fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructStencilDestroy,(ex->hstencil));
180f4bdf6c4SBarry Smith   if (ex->needsinitialization) {
181fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_StructMatrixInitialize,(ex->hmat));
182f4bdf6c4SBarry Smith     ex->needsinitialization = PETSC_FALSE;
183f4bdf6c4SBarry Smith   }
184f4bdf6c4SBarry Smith 
185f4bdf6c4SBarry Smith   /* set the global and local sizes of the matrix */
186f4bdf6c4SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
187f4bdf6c4SBarry Smith   ierr = MatSetSizes(mat,dof*nx*ny*nz,dof*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
188f4bdf6c4SBarry Smith   ierr = PetscLayoutSetBlockSize(mat->rmap,1);CHKERRQ(ierr);
189f4bdf6c4SBarry Smith   ierr = PetscLayoutSetBlockSize(mat->cmap,1);CHKERRQ(ierr);
190f4bdf6c4SBarry Smith   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
191f4bdf6c4SBarry Smith   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
192f4bdf6c4SBarry Smith 
193f4bdf6c4SBarry Smith   if (dim == 3) {
194f4bdf6c4SBarry Smith     mat->ops->setvalueslocal = MatSetValuesLocal_HYPREStruct_3d;
195f4bdf6c4SBarry Smith     mat->ops->zerorowslocal  = MatZeroRowsLocal_HYPREStruct_3d;
196f4bdf6c4SBarry Smith     mat->ops->zeroentries    = MatZeroEntries_HYPREStruct_3d;
1978865f1eaSKarl Rupp 
198f4bdf6c4SBarry Smith     ierr = MatZeroEntries_HYPREStruct_3d(mat);CHKERRQ(ierr);
199ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only support for 3d DMDA currently");
200f4bdf6c4SBarry Smith 
201f4bdf6c4SBarry Smith   /* get values that will be used repeatedly in MatSetValuesLocal() and MatZeroRowsLocal() repeatedly */
2020298fd71SBarry Smith   ierr        = MatGetOwnershipRange(mat,&ex->rstart,NULL);CHKERRQ(ierr);
203565245c5SBarry Smith   ierr        = DMGetLocalToGlobalMapping(ex->da,&ltog);CHKERRQ(ierr);
204565245c5SBarry Smith   ierr        = ISLocalToGlobalMappingGetIndices(ltog, (const PetscInt **) &ex->gindices);CHKERRQ(ierr);
205f4bdf6c4SBarry Smith   ierr        = DMDAGetGhostCorners(ex->da,0,0,0,&ex->gnx,&ex->gnxgny,0);CHKERRQ(ierr);
206f4bdf6c4SBarry Smith   ex->gnxgny *= ex->gnx;
207f4bdf6c4SBarry Smith   ierr        = DMDAGetCorners(ex->da,&ex->xs,&ex->ys,&ex->zs,&ex->nx,&ex->ny,0);CHKERRQ(ierr);
208f4bdf6c4SBarry Smith   ex->nxny    = ex->nx*ex->ny;
209f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
210f4bdf6c4SBarry Smith }
211f4bdf6c4SBarry Smith 
212f4bdf6c4SBarry Smith PetscErrorCode MatMult_HYPREStruct(Mat A,Vec x,Vec y)
213f4bdf6c4SBarry Smith {
214f4bdf6c4SBarry Smith   PetscErrorCode    ierr;
215*b026d285SBarry Smith   const PetscScalar *xx;
216*b026d285SBarry Smith   PetscScalar       *yy;
2174ddd07fcSJed Brown   PetscInt          ilower[3],iupper[3];
218f4bdf6c4SBarry Smith   Mat_HYPREStruct   *mx = (Mat_HYPREStruct*)(A->data);
219f4bdf6c4SBarry Smith 
220f4bdf6c4SBarry Smith   PetscFunctionBegin;
221f4bdf6c4SBarry Smith   ierr = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
2228865f1eaSKarl Rupp 
223f4bdf6c4SBarry Smith   iupper[0] += ilower[0] - 1;
224f4bdf6c4SBarry Smith   iupper[1] += ilower[1] - 1;
225f4bdf6c4SBarry Smith   iupper[2] += ilower[2] - 1;
226f4bdf6c4SBarry Smith 
227f4bdf6c4SBarry Smith   /* copy x values over to hypre */
228fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorSetConstantValues,(mx->hb,0.0));
229*b026d285SBarry Smith   ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
230*b026d285SBarry Smith   PetscStackCallStandard(HYPRE_StructVectorSetBoxValues,(mx->hb,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,(PetscScalar*)xx));
231*b026d285SBarry Smith   ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
232fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorAssemble,(mx->hb));
233fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixMatvec,(1.0,mx->hmat,mx->hb,0.0,mx->hx));
234f4bdf6c4SBarry Smith 
235f4bdf6c4SBarry Smith   /* copy solution values back to PETSc */
236f4bdf6c4SBarry Smith   ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
2378b1f7689SBarry Smith   PetscStackCallStandard(HYPRE_StructVectorGetBoxValues,(mx->hx,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,yy));
238f4bdf6c4SBarry Smith   ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
239f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
240f4bdf6c4SBarry Smith }
241f4bdf6c4SBarry Smith 
242f4bdf6c4SBarry Smith PetscErrorCode MatAssemblyEnd_HYPREStruct(Mat mat,MatAssemblyType mode)
243f4bdf6c4SBarry Smith {
244f4bdf6c4SBarry Smith   Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data;
245f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
246f4bdf6c4SBarry Smith 
247f4bdf6c4SBarry Smith   PetscFunctionBegin;
248fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixAssemble,(ex->hmat));
249fd3f9acdSBarry Smith   /* PetscStackCallStandard(HYPRE_StructMatrixPrint,("dummy",ex->hmat,0)); */
250f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
251f4bdf6c4SBarry Smith }
252f4bdf6c4SBarry Smith 
253f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPREStruct(Mat mat)
254f4bdf6c4SBarry Smith {
255f4bdf6c4SBarry Smith   PetscFunctionBegin;
256f4bdf6c4SBarry Smith   /* before the DMDA is set to the matrix the zero doesn't need to do anything */
257f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
258f4bdf6c4SBarry Smith }
259f4bdf6c4SBarry Smith 
260f4bdf6c4SBarry Smith PetscErrorCode MatDestroy_HYPREStruct(Mat mat)
261f4bdf6c4SBarry Smith {
262f4bdf6c4SBarry Smith   Mat_HYPREStruct *ex = (Mat_HYPREStruct*) mat->data;
263f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
264f4bdf6c4SBarry Smith 
265f4bdf6c4SBarry Smith   PetscFunctionBegin;
266fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructMatrixDestroy,(ex->hmat));
267fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorDestroy,(ex->hx));
268fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_StructVectorDestroy,(ex->hb));
269f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
270f4bdf6c4SBarry Smith }
271f4bdf6c4SBarry Smith 
2728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_HYPREStruct(Mat B)
273f4bdf6c4SBarry Smith {
274f4bdf6c4SBarry Smith   Mat_HYPREStruct *ex;
275f4bdf6c4SBarry Smith   PetscErrorCode  ierr;
276f4bdf6c4SBarry Smith 
277f4bdf6c4SBarry Smith   PetscFunctionBegin;
278b00a9115SJed Brown   ierr         = PetscNewLog(B,&ex);CHKERRQ(ierr);
279f4bdf6c4SBarry Smith   B->data      = (void*)ex;
280f4bdf6c4SBarry Smith   B->rmap->bs  = 1;
281f4bdf6c4SBarry Smith   B->assembled = PETSC_FALSE;
282f4bdf6c4SBarry Smith 
283f4bdf6c4SBarry Smith   B->insertmode = NOT_SET_VALUES;
284f4bdf6c4SBarry Smith 
285f4bdf6c4SBarry Smith   B->ops->assemblyend = MatAssemblyEnd_HYPREStruct;
286f4bdf6c4SBarry Smith   B->ops->mult        = MatMult_HYPREStruct;
287f4bdf6c4SBarry Smith   B->ops->zeroentries = MatZeroEntries_HYPREStruct;
288f4bdf6c4SBarry Smith   B->ops->destroy     = MatDestroy_HYPREStruct;
289f4bdf6c4SBarry Smith 
290f4bdf6c4SBarry Smith   ex->needsinitialization = PETSC_TRUE;
291f4bdf6c4SBarry Smith 
292ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&(ex->hcomm));CHKERRQ(ierr);
293bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSetupDM_C",MatSetupDM_HYPREStruct);CHKERRQ(ierr);
294f4bdf6c4SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRESTRUCT);CHKERRQ(ierr);
295f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
296f4bdf6c4SBarry Smith }
297f4bdf6c4SBarry Smith 
298f4bdf6c4SBarry Smith /*MC
299f4bdf6c4SBarry Smith    MATHYPRESSTRUCT - MATHYPRESSTRUCT = "hypresstruct" - A matrix type to be used for parallel sparse matrices
300f4bdf6c4SBarry Smith           based on the hypre HYPRE_SStructMatrix.
301f4bdf6c4SBarry Smith 
302f4bdf6c4SBarry Smith 
303f4bdf6c4SBarry Smith    Level: intermediate
304f4bdf6c4SBarry Smith 
305f4bdf6c4SBarry Smith    Notes: Unlike hypre's general semi-struct object consisting of a collection of structured-grid objects and unstructured
306*b026d285SBarry Smith           grid objects, we restrict the semi-struct objects to consist of only structured-grid components.
307f4bdf6c4SBarry Smith 
308f4bdf6c4SBarry 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
309f4bdf6c4SBarry Smith           be defined by a DMDA.
310f4bdf6c4SBarry Smith 
311c688c046SMatthew G Knepley           The matrix needs a DMDA associated with it by either a call to MatSetupDM() or if the matrix is obtained from DMCreateMatrix()
312f4bdf6c4SBarry Smith 
313f4bdf6c4SBarry Smith M*/
314f4bdf6c4SBarry Smith 
315f4bdf6c4SBarry Smith PetscErrorCode  MatSetValuesLocal_HYPRESStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscInt ncol,const PetscInt icol[],const PetscScalar y[],InsertMode addv)
316f4bdf6c4SBarry Smith {
317f4bdf6c4SBarry Smith   PetscErrorCode    ierr;
318f4bdf6c4SBarry Smith   PetscInt          i,j,stencil,index[3];
319f4bdf6c4SBarry Smith   const PetscScalar *values = y;
320f4bdf6c4SBarry Smith   Mat_HYPRESStruct  *ex     = (Mat_HYPRESStruct*) mat->data;
321f4bdf6c4SBarry Smith 
3224ddd07fcSJed Brown   PetscInt          part = 0;          /* Petsc sstruct interface only allows 1 part */
3234ddd07fcSJed Brown   PetscInt          ordering;
3244ddd07fcSJed Brown   PetscInt          grid_rank, to_grid_rank;
3254ddd07fcSJed Brown   PetscInt          var_type, to_var_type;
3264ddd07fcSJed Brown   PetscInt          to_var_entry = 0;
3274ddd07fcSJed Brown   PetscInt          nvars= ex->nvars;
328f4bdf6c4SBarry Smith   PetscInt          row,*entries;
329f4bdf6c4SBarry Smith 
330f4bdf6c4SBarry Smith   PetscFunctionBegin;
331785e854fSJed Brown   ierr     = PetscMalloc1(7*nvars,&entries);CHKERRQ(ierr);
332f4bdf6c4SBarry Smith   ordering = ex->dofs_order;  /* ordering= 0   nodal ordering
333f4bdf6c4SBarry Smith                                            1   variable ordering */
334f4bdf6c4SBarry Smith   /* stencil entries are orderer by variables: var0_stencil0, var0_stencil1, ..., var0_stencil6, var1_stencil0, var1_stencil1, ...  */
335f4bdf6c4SBarry Smith 
336f4bdf6c4SBarry Smith   if (!ordering) {
337f4bdf6c4SBarry Smith     for (i=0; i<nrow; i++) {
338f4bdf6c4SBarry Smith       grid_rank= irow[i]/nvars;
339f4bdf6c4SBarry Smith       var_type = (irow[i] % nvars);
340f4bdf6c4SBarry Smith 
341f4bdf6c4SBarry Smith       for (j=0; j<ncol; j++) {
342f4bdf6c4SBarry Smith         to_grid_rank= icol[j]/nvars;
343f4bdf6c4SBarry Smith         to_var_type = (icol[j] % nvars);
344f4bdf6c4SBarry Smith 
345f4bdf6c4SBarry Smith         to_var_entry= to_var_entry*7;
346f4bdf6c4SBarry Smith         entries[j]  = to_var_entry;
347f4bdf6c4SBarry Smith 
348f4bdf6c4SBarry Smith         stencil = to_grid_rank-grid_rank;
349f4bdf6c4SBarry Smith         if (!stencil) {
350f4bdf6c4SBarry Smith           entries[j] += 3;
351f4bdf6c4SBarry Smith         } else if (stencil == -1) {
352f4bdf6c4SBarry Smith           entries[j] += 2;
353f4bdf6c4SBarry Smith         } else if (stencil == 1) {
354f4bdf6c4SBarry Smith           entries[j] += 4;
355f4bdf6c4SBarry Smith         } else if (stencil == -ex->gnx) {
356f4bdf6c4SBarry Smith           entries[j] += 1;
357f4bdf6c4SBarry Smith         } else if (stencil == ex->gnx) {
358f4bdf6c4SBarry Smith           entries[j] += 5;
359f4bdf6c4SBarry Smith         } else if (stencil == -ex->gnxgny) {
360f4bdf6c4SBarry Smith           entries[j] += 0;
361f4bdf6c4SBarry Smith         } else if (stencil == ex->gnxgny) {
362f4bdf6c4SBarry Smith           entries[j] += 6;
363f4bdf6c4SBarry 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);
364f4bdf6c4SBarry Smith       }
365f4bdf6c4SBarry Smith 
366f4bdf6c4SBarry Smith       row      = ex->gindices[grid_rank] - ex->rstart;
367f4bdf6c4SBarry Smith       index[0] = ex->xs + (row % ex->nx);
368f4bdf6c4SBarry Smith       index[1] = ex->ys + ((row/ex->nx) % ex->ny);
369f4bdf6c4SBarry Smith       index[2] = ex->zs + (row/(ex->nxny));
370f4bdf6c4SBarry Smith 
371f4bdf6c4SBarry Smith       if (addv == ADD_VALUES) {
3728b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructMatrixAddToValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
373f4bdf6c4SBarry Smith       } else {
3748b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
375f4bdf6c4SBarry Smith       }
376f4bdf6c4SBarry Smith       values += ncol;
377f4bdf6c4SBarry Smith     }
378f4bdf6c4SBarry Smith   } else {
379f4bdf6c4SBarry Smith     for (i=0; i<nrow; i++) {
380f4bdf6c4SBarry Smith       var_type = irow[i]/(ex->gnxgnygnz);
381f4bdf6c4SBarry Smith       grid_rank= irow[i] - var_type*(ex->gnxgnygnz);
382f4bdf6c4SBarry Smith 
383f4bdf6c4SBarry Smith       for (j=0; j<ncol; j++) {
384f4bdf6c4SBarry Smith         to_var_type = icol[j]/(ex->gnxgnygnz);
385f4bdf6c4SBarry Smith         to_grid_rank= icol[j] - to_var_type*(ex->gnxgnygnz);
386f4bdf6c4SBarry Smith 
387f4bdf6c4SBarry Smith         to_var_entry= to_var_entry*7;
388f4bdf6c4SBarry Smith         entries[j]  = to_var_entry;
389f4bdf6c4SBarry Smith 
390f4bdf6c4SBarry Smith         stencil = to_grid_rank-grid_rank;
391f4bdf6c4SBarry Smith         if (!stencil) {
392f4bdf6c4SBarry Smith           entries[j] += 3;
393f4bdf6c4SBarry Smith         } else if (stencil == -1) {
394f4bdf6c4SBarry Smith           entries[j] += 2;
395f4bdf6c4SBarry Smith         } else if (stencil == 1) {
396f4bdf6c4SBarry Smith           entries[j] += 4;
397f4bdf6c4SBarry Smith         } else if (stencil == -ex->gnx) {
398f4bdf6c4SBarry Smith           entries[j] += 1;
399f4bdf6c4SBarry Smith         } else if (stencil == ex->gnx) {
400f4bdf6c4SBarry Smith           entries[j] += 5;
401f4bdf6c4SBarry Smith         } else if (stencil == -ex->gnxgny) {
402f4bdf6c4SBarry Smith           entries[j] += 0;
403f4bdf6c4SBarry Smith         } else if (stencil == ex->gnxgny) {
404f4bdf6c4SBarry Smith           entries[j] += 6;
405f4bdf6c4SBarry 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);
406f4bdf6c4SBarry Smith       }
407f4bdf6c4SBarry Smith 
408f4bdf6c4SBarry Smith       row      = ex->gindices[grid_rank] - ex->rstart;
409f4bdf6c4SBarry Smith       index[0] = ex->xs + (row % ex->nx);
410f4bdf6c4SBarry Smith       index[1] = ex->ys + ((row/ex->nx) % ex->ny);
411f4bdf6c4SBarry Smith       index[2] = ex->zs + (row/(ex->nxny));
412f4bdf6c4SBarry Smith 
413f4bdf6c4SBarry Smith       if (addv == ADD_VALUES) {
4148b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructMatrixAddToValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
415f4bdf6c4SBarry Smith       } else {
4168b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,ncol,(HYPRE_Int *)entries,(PetscScalar*)values));
417f4bdf6c4SBarry Smith       }
418f4bdf6c4SBarry Smith       values += ncol;
419f4bdf6c4SBarry Smith     }
420f4bdf6c4SBarry Smith   }
421f4bdf6c4SBarry Smith   ierr = PetscFree(entries);CHKERRQ(ierr);
422f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
423f4bdf6c4SBarry Smith }
424f4bdf6c4SBarry Smith 
425f4bdf6c4SBarry Smith PetscErrorCode  MatZeroRowsLocal_HYPRESStruct_3d(Mat mat,PetscInt nrow,const PetscInt irow[],PetscScalar d,Vec x,Vec b)
426f4bdf6c4SBarry Smith {
427f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
428f4bdf6c4SBarry Smith   PetscInt         i,index[3];
429f4bdf6c4SBarry Smith   PetscScalar      **values;
430f4bdf6c4SBarry Smith   Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data;
431f4bdf6c4SBarry Smith 
4324ddd07fcSJed Brown   PetscInt         part     = 0;      /* Petsc sstruct interface only allows 1 part */
4334ddd07fcSJed Brown   PetscInt         ordering = ex->dofs_order;
4344ddd07fcSJed Brown   PetscInt         grid_rank;
4354ddd07fcSJed Brown   PetscInt         var_type;
4364ddd07fcSJed Brown   PetscInt         nvars= ex->nvars;
437f4bdf6c4SBarry Smith   PetscInt         row,*entries;
438f4bdf6c4SBarry Smith 
439f4bdf6c4SBarry Smith   PetscFunctionBegin;
440ce94432eSBarry Smith   if (x && b) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"No support");
441785e854fSJed Brown   ierr = PetscMalloc1(7*nvars,&entries);CHKERRQ(ierr);
442f4bdf6c4SBarry Smith 
443785e854fSJed Brown   ierr = PetscMalloc1(nvars,&values);CHKERRQ(ierr);
444785e854fSJed Brown   ierr = PetscMalloc1(7*nvars*nvars,&values[0]);CHKERRQ(ierr);
445f4bdf6c4SBarry Smith   for (i=1; i<nvars; i++) {
446f4bdf6c4SBarry Smith     values[i] = values[i-1] + nvars*7;
447f4bdf6c4SBarry Smith   }
448f4bdf6c4SBarry Smith 
449f4bdf6c4SBarry Smith   for (i=0; i< nvars; i++) {
450f4bdf6c4SBarry Smith     ierr           = PetscMemzero(values[i],nvars*7*sizeof(PetscScalar));CHKERRQ(ierr);
451f4bdf6c4SBarry Smith     *(values[i]+3) = d;
452f4bdf6c4SBarry Smith   }
453f4bdf6c4SBarry Smith 
4548865f1eaSKarl Rupp   for (i= 0; i< nvars*7; i++) entries[i] = i;
455f4bdf6c4SBarry Smith 
456f4bdf6c4SBarry Smith   if (!ordering) {
457f4bdf6c4SBarry Smith     for (i=0; i<nrow; i++) {
458f4bdf6c4SBarry Smith       grid_rank = irow[i] / nvars;
459f4bdf6c4SBarry Smith       var_type  =(irow[i] % nvars);
460f4bdf6c4SBarry Smith 
461f4bdf6c4SBarry Smith       row      = ex->gindices[grid_rank] - ex->rstart;
462f4bdf6c4SBarry Smith       index[0] = ex->xs + (row % ex->nx);
463f4bdf6c4SBarry Smith       index[1] = ex->ys + ((row/ex->nx) % ex->ny);
464f4bdf6c4SBarry Smith       index[2] = ex->zs + (row/(ex->nxny));
4658b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,7*nvars,(HYPRE_Int *)entries,values[var_type]));
466f4bdf6c4SBarry Smith     }
467f4bdf6c4SBarry Smith   } else {
468f4bdf6c4SBarry Smith     for (i=0; i<nrow; i++) {
469f4bdf6c4SBarry Smith       var_type = irow[i]/(ex->gnxgnygnz);
470f4bdf6c4SBarry Smith       grid_rank= irow[i] - var_type*(ex->gnxgnygnz);
471f4bdf6c4SBarry Smith 
472f4bdf6c4SBarry Smith       row      = ex->gindices[grid_rank] - ex->rstart;
473f4bdf6c4SBarry Smith       index[0] = ex->xs + (row % ex->nx);
474f4bdf6c4SBarry Smith       index[1] = ex->ys + ((row/ex->nx) % ex->ny);
475f4bdf6c4SBarry Smith       index[2] = ex->zs + (row/(ex->nxny));
4768b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructMatrixSetValues,(ex->ss_mat,part,(HYPRE_Int *)index,var_type,7*nvars,(HYPRE_Int *)entries,values[var_type]));
477f4bdf6c4SBarry Smith     }
478f4bdf6c4SBarry Smith   }
479fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructMatrixAssemble,(ex->ss_mat));
480f4bdf6c4SBarry Smith   ierr = PetscFree(values[0]);CHKERRQ(ierr);
481f4bdf6c4SBarry Smith   ierr = PetscFree(values);CHKERRQ(ierr);
482f4bdf6c4SBarry Smith   ierr = PetscFree(entries);CHKERRQ(ierr);
483f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
484f4bdf6c4SBarry Smith }
485f4bdf6c4SBarry Smith 
486f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPRESStruct_3d(Mat mat)
487f4bdf6c4SBarry Smith {
488f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
489f4bdf6c4SBarry Smith   Mat_HYPRESStruct *ex  = (Mat_HYPRESStruct*) mat->data;
4904ddd07fcSJed Brown   PetscInt         nvars= ex->nvars;
4914ddd07fcSJed Brown   PetscInt         size;
4924ddd07fcSJed Brown   PetscInt         part= 0;   /* only one part */
493f4bdf6c4SBarry Smith 
494f4bdf6c4SBarry Smith   PetscFunctionBegin;
495f4bdf6c4SBarry 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);
496f4bdf6c4SBarry Smith   {
4974ddd07fcSJed Brown     PetscInt    i,*entries,iupper[3],ilower[3];
498f4bdf6c4SBarry Smith     PetscScalar *values;
4994ddd07fcSJed Brown 
500f4bdf6c4SBarry Smith 
501f4bdf6c4SBarry Smith     for (i= 0; i< 3; i++) {
502f4bdf6c4SBarry Smith       ilower[i]= ex->hbox.imin[i];
503f4bdf6c4SBarry Smith       iupper[i]= ex->hbox.imax[i];
504f4bdf6c4SBarry Smith     }
505f4bdf6c4SBarry Smith 
506dcca6d9dSJed Brown     ierr = PetscMalloc2(nvars*7,&entries,nvars*7*size,&values);CHKERRQ(ierr);
5078865f1eaSKarl Rupp     for (i= 0; i< nvars*7; i++) entries[i]= i;
508f4bdf6c4SBarry Smith     ierr = PetscMemzero(values,nvars*7*size*sizeof(PetscScalar));CHKERRQ(ierr);
509f4bdf6c4SBarry Smith 
510f4bdf6c4SBarry Smith     for (i= 0; i< nvars; i++) {
5118b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructMatrixSetBoxValues,(ex->ss_mat,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,nvars*7,(HYPRE_Int *)entries,values));
512f4bdf6c4SBarry Smith     }
513f4bdf6c4SBarry Smith     ierr = PetscFree2(entries,values);CHKERRQ(ierr);
514f4bdf6c4SBarry Smith   }
515fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructMatrixAssemble,(ex->ss_mat));
516f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
517f4bdf6c4SBarry Smith }
518f4bdf6c4SBarry Smith 
519f4bdf6c4SBarry Smith 
5201e6b0712SBarry Smith static PetscErrorCode  MatSetupDM_HYPRESStruct(Mat mat,DM da)
521f4bdf6c4SBarry Smith {
522f4bdf6c4SBarry Smith   PetscErrorCode         ierr;
523f4bdf6c4SBarry Smith   Mat_HYPRESStruct       *ex = (Mat_HYPRESStruct*) mat->data;
524f4bdf6c4SBarry Smith   PetscInt               dim,dof,sw[3],nx,ny,nz;
5254ddd07fcSJed Brown   PetscInt               ilower[3],iupper[3],ssize,i;
526bff4a2f0SMatthew G. Knepley   DMBoundaryType         px,py,pz;
527f4bdf6c4SBarry Smith   DMDAStencilType        st;
5284ddd07fcSJed Brown   PetscInt               nparts= 1;  /* assuming only one part */
5294ddd07fcSJed Brown   PetscInt               part  = 0;
530565245c5SBarry Smith   ISLocalToGlobalMapping ltog;
531*b026d285SBarry Smith 
532f4bdf6c4SBarry Smith   PetscFunctionBegin;
533f4bdf6c4SBarry Smith   ex->da = da;
534f4bdf6c4SBarry Smith   ierr   = PetscObjectReference((PetscObject)da);CHKERRQ(ierr);
535f4bdf6c4SBarry Smith 
5367ab44359SJed Brown   ierr       = DMDAGetInfo(ex->da,&dim,0,0,0,0,0,0,&dof,&sw[0],&px,&py,&pz,&st);CHKERRQ(ierr);
537f4bdf6c4SBarry Smith   ierr       = DMDAGetCorners(ex->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
538f4bdf6c4SBarry Smith   iupper[0] += ilower[0] - 1;
539f4bdf6c4SBarry Smith   iupper[1] += ilower[1] - 1;
540f4bdf6c4SBarry Smith   iupper[2] += ilower[2] - 1;
541f4bdf6c4SBarry Smith   /* the hypre_Box is used to zero out the matrix entries in MatZeroValues() */
542f4bdf6c4SBarry Smith   ex->hbox.imin[0] = ilower[0];
543f4bdf6c4SBarry Smith   ex->hbox.imin[1] = ilower[1];
544f4bdf6c4SBarry Smith   ex->hbox.imin[2] = ilower[2];
545f4bdf6c4SBarry Smith   ex->hbox.imax[0] = iupper[0];
546f4bdf6c4SBarry Smith   ex->hbox.imax[1] = iupper[1];
547f4bdf6c4SBarry Smith   ex->hbox.imax[2] = iupper[2];
548f4bdf6c4SBarry Smith 
549f4bdf6c4SBarry Smith   ex->dofs_order = 0;
550f4bdf6c4SBarry Smith 
551f4bdf6c4SBarry Smith   /* assuming that the same number of dofs on each gridpoint. Also assume all cell-centred based */
552f4bdf6c4SBarry Smith   ex->nvars= dof;
553f4bdf6c4SBarry Smith 
554f4bdf6c4SBarry Smith   /* create the hypre grid object and set its information */
555ce94432eSBarry Smith   if (px || py || pz) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add periodic support by calling HYPRE_SStructGridSetPeriodic()");
556fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGridCreate,(ex->hcomm,dim,nparts,&ex->ss_grid));
557fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGridSetExtents,(ex->ss_grid,part,ex->hbox.imin,ex->hbox.imax));
558f4bdf6c4SBarry Smith   {
559f4bdf6c4SBarry Smith     HYPRE_SStructVariable *vartypes;
560785e854fSJed Brown     ierr = PetscMalloc1(ex->nvars,&vartypes);CHKERRQ(ierr);
5618865f1eaSKarl Rupp     for (i= 0; i< ex->nvars; i++) vartypes[i]= HYPRE_SSTRUCT_VARIABLE_CELL;
562fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructGridSetVariables,(ex->ss_grid, part, ex->nvars,vartypes));
563f4bdf6c4SBarry Smith     ierr = PetscFree(vartypes);CHKERRQ(ierr);
564f4bdf6c4SBarry Smith   }
565fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGridAssemble,(ex->ss_grid));
566f4bdf6c4SBarry Smith 
567f4bdf6c4SBarry Smith   sw[1] = sw[0];
568f4bdf6c4SBarry Smith   sw[2] = sw[1];
569fd3f9acdSBarry Smith   /* PetscStackCallStandard(HYPRE_SStructGridSetNumGhost,(ex->ss_grid,sw)); */
570f4bdf6c4SBarry Smith 
571f4bdf6c4SBarry Smith   /* create the hypre stencil object and set its information */
572ce94432eSBarry Smith   if (sw[0] > 1) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add support for wider stencils");
573ce94432eSBarry Smith   if (st == DMDA_STENCIL_BOX) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Ask us to add support for box stencils");
574f4bdf6c4SBarry Smith 
575f4bdf6c4SBarry Smith   if (dim == 1) {
5764ddd07fcSJed Brown     PetscInt offsets[3][1] = {{-1},{0},{1}};
5774ddd07fcSJed Brown     PetscInt j, cnt;
578f4bdf6c4SBarry Smith 
579f4bdf6c4SBarry Smith     ssize = 3*(ex->nvars);
580fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil));
581f4bdf6c4SBarry Smith     cnt= 0;
582f4bdf6c4SBarry Smith     for (i = 0; i < (ex->nvars); i++) {
583f4bdf6c4SBarry Smith       for (j = 0; j < 3; j++) {
5848b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, (HYPRE_Int *)offsets[j], i));
585f4bdf6c4SBarry Smith         cnt++;
586f4bdf6c4SBarry Smith       }
587f4bdf6c4SBarry Smith     }
588f4bdf6c4SBarry Smith   } else if (dim == 2) {
5894ddd07fcSJed Brown     PetscInt offsets[5][2] = {{0,-1},{-1,0},{0,0},{1,0},{0,1}};
5904ddd07fcSJed Brown     PetscInt j, cnt;
591f4bdf6c4SBarry Smith 
592f4bdf6c4SBarry Smith     ssize = 5*(ex->nvars);
593fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil));
594f4bdf6c4SBarry Smith     cnt= 0;
595f4bdf6c4SBarry Smith     for (i= 0; i< (ex->nvars); i++) {
596f4bdf6c4SBarry Smith       for (j= 0; j< 5; j++) {
5978b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, (HYPRE_Int *)offsets[j], i));
598f4bdf6c4SBarry Smith         cnt++;
599f4bdf6c4SBarry Smith       }
600f4bdf6c4SBarry Smith     }
601f4bdf6c4SBarry Smith   } else if (dim == 3) {
6024ddd07fcSJed 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}};
6034ddd07fcSJed Brown     PetscInt j, cnt;
604f4bdf6c4SBarry Smith 
605f4bdf6c4SBarry Smith     ssize = 7*(ex->nvars);
606fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructStencilCreate,(dim,ssize,&ex->ss_stencil));
607f4bdf6c4SBarry Smith     cnt= 0;
608f4bdf6c4SBarry Smith     for (i= 0; i< (ex->nvars); i++) {
609f4bdf6c4SBarry Smith       for (j= 0; j< 7; j++) {
6108b1f7689SBarry Smith         PetscStackCallStandard(HYPRE_SStructStencilSetEntry,(ex->ss_stencil, cnt, (HYPRE_Int *)offsets[j], i));
611f4bdf6c4SBarry Smith         cnt++;
612f4bdf6c4SBarry Smith       }
613f4bdf6c4SBarry Smith     }
614f4bdf6c4SBarry Smith   }
615f4bdf6c4SBarry Smith 
616f4bdf6c4SBarry Smith   /* create the HYPRE graph */
617fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGraphCreate,(ex->hcomm, ex->ss_grid, &(ex->ss_graph)));
618f4bdf6c4SBarry Smith 
619f4bdf6c4SBarry Smith   /* set the stencil graph. Note that each variable has the same graph. This means that each
620f4bdf6c4SBarry Smith      variable couples to all the other variable and with the same stencil pattern. */
621f4bdf6c4SBarry Smith   for (i= 0; i< (ex->nvars); i++) {
622fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructGraphSetStencil,(ex->ss_graph,part,i,ex->ss_stencil));
623f4bdf6c4SBarry Smith   }
624fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGraphAssemble,(ex->ss_graph));
625f4bdf6c4SBarry Smith 
626f4bdf6c4SBarry Smith   /* create the HYPRE sstruct vectors for rhs and solution */
627fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorCreate,(ex->hcomm,ex->ss_grid,&ex->ss_b));
628fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorCreate,(ex->hcomm,ex->ss_grid,&ex->ss_x));
629fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorInitialize,(ex->ss_b));
630fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorInitialize,(ex->ss_x));
631fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorAssemble,(ex->ss_b));
632fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorAssemble,(ex->ss_x));
633f4bdf6c4SBarry Smith 
634f4bdf6c4SBarry Smith   /* create the hypre matrix object and set its information */
635fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructMatrixCreate,(ex->hcomm,ex->ss_graph,&ex->ss_mat));
636fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGridDestroy,(ex->ss_grid));
637fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructStencilDestroy,(ex->ss_stencil));
638f4bdf6c4SBarry Smith   if (ex->needsinitialization) {
639fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructMatrixInitialize,(ex->ss_mat));
640f4bdf6c4SBarry Smith     ex->needsinitialization = PETSC_FALSE;
641f4bdf6c4SBarry Smith   }
642f4bdf6c4SBarry Smith 
643f4bdf6c4SBarry Smith   /* set the global and local sizes of the matrix */
644f4bdf6c4SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
645f4bdf6c4SBarry Smith   ierr = MatSetSizes(mat,dof*nx*ny*nz,dof*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
646f4bdf6c4SBarry Smith   ierr = PetscLayoutSetBlockSize(mat->rmap,1);CHKERRQ(ierr);
647f4bdf6c4SBarry Smith   ierr = PetscLayoutSetBlockSize(mat->cmap,1);CHKERRQ(ierr);
648f4bdf6c4SBarry Smith   ierr = PetscLayoutSetUp(mat->rmap);CHKERRQ(ierr);
649f4bdf6c4SBarry Smith   ierr = PetscLayoutSetUp(mat->cmap);CHKERRQ(ierr);
650f4bdf6c4SBarry Smith 
651f4bdf6c4SBarry Smith   if (dim == 3) {
652f4bdf6c4SBarry Smith     mat->ops->setvalueslocal = MatSetValuesLocal_HYPRESStruct_3d;
653f4bdf6c4SBarry Smith     mat->ops->zerorowslocal  = MatZeroRowsLocal_HYPRESStruct_3d;
654f4bdf6c4SBarry Smith     mat->ops->zeroentries    = MatZeroEntries_HYPRESStruct_3d;
6558865f1eaSKarl Rupp 
656f4bdf6c4SBarry Smith     ierr = MatZeroEntries_HYPRESStruct_3d(mat);CHKERRQ(ierr);
657ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Only support for 3d DMDA currently");
658f4bdf6c4SBarry Smith 
659f4bdf6c4SBarry Smith   /* get values that will be used repeatedly in MatSetValuesLocal() and MatZeroRowsLocal() repeatedly */
6600298fd71SBarry Smith   ierr = MatGetOwnershipRange(mat,&ex->rstart,NULL);CHKERRQ(ierr);
661565245c5SBarry Smith   ierr = DMGetLocalToGlobalMapping(ex->da,&ltog);CHKERRQ(ierr);
662565245c5SBarry Smith   ierr = ISLocalToGlobalMappingGetIndices(ltog, (const PetscInt **) &ex->gindices);CHKERRQ(ierr);
663f4bdf6c4SBarry Smith   ierr = DMDAGetGhostCorners(ex->da,0,0,0,&ex->gnx,&ex->gnxgny,&ex->gnxgnygnz);CHKERRQ(ierr);
6648865f1eaSKarl Rupp 
665f4bdf6c4SBarry Smith   ex->gnxgny    *= ex->gnx;
666f4bdf6c4SBarry Smith   ex->gnxgnygnz *= ex->gnxgny;
6678865f1eaSKarl Rupp 
668f4bdf6c4SBarry Smith   ierr = DMDAGetCorners(ex->da,&ex->xs,&ex->ys,&ex->zs,&ex->nx,&ex->ny,&ex->nz);CHKERRQ(ierr);
6698865f1eaSKarl Rupp 
670f4bdf6c4SBarry Smith   ex->nxny   = ex->nx*ex->ny;
671f4bdf6c4SBarry Smith   ex->nxnynz = ex->nz*ex->nxny;
672f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
673f4bdf6c4SBarry Smith }
674f4bdf6c4SBarry Smith 
675f4bdf6c4SBarry Smith PetscErrorCode MatMult_HYPRESStruct(Mat A,Vec x,Vec y)
676f4bdf6c4SBarry Smith {
677f4bdf6c4SBarry Smith   PetscErrorCode    ierr;
678*b026d285SBarry Smith   const PetscScalar *xx;
679*b026d285SBarry Smith   PetscScalar       *yy;
6804ddd07fcSJed Brown   PetscInt          ilower[3],iupper[3];
681f4bdf6c4SBarry Smith   Mat_HYPRESStruct  *mx     = (Mat_HYPRESStruct*)(A->data);
6824ddd07fcSJed Brown   PetscInt          ordering= mx->dofs_order;
6834ddd07fcSJed Brown   PetscInt          nvars   = mx->nvars;
6844ddd07fcSJed Brown   PetscInt          part    = 0;
6854ddd07fcSJed Brown   PetscInt          size;
6864ddd07fcSJed Brown   PetscInt          i;
687f4bdf6c4SBarry Smith 
688f4bdf6c4SBarry Smith   PetscFunctionBegin;
689f4bdf6c4SBarry Smith   ierr       = DMDAGetCorners(mx->da,&ilower[0],&ilower[1],&ilower[2],&iupper[0],&iupper[1],&iupper[2]);CHKERRQ(ierr);
690f4bdf6c4SBarry Smith   iupper[0] += ilower[0] - 1;
691f4bdf6c4SBarry Smith   iupper[1] += ilower[1] - 1;
692f4bdf6c4SBarry Smith   iupper[2] += ilower[2] - 1;
693f4bdf6c4SBarry Smith 
694f4bdf6c4SBarry Smith   size= 1;
6958865f1eaSKarl Rupp   for (i = 0; i < 3; i++) size *= (iupper[i]-ilower[i]+1);
696f4bdf6c4SBarry Smith 
697f4bdf6c4SBarry Smith   /* copy x values over to hypre for variable ordering */
698f4bdf6c4SBarry Smith   if (ordering) {
699fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
700*b026d285SBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
701f4bdf6c4SBarry Smith     for (i= 0; i< nvars; i++) {
702*b026d285SBarry Smith       PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,(PetscScalar*)xx+(size*i)));
703f4bdf6c4SBarry Smith     }
704*b026d285SBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
705fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
706fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
707f4bdf6c4SBarry Smith 
708f4bdf6c4SBarry Smith     /* copy solution values back to PETSc */
709f4bdf6c4SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
710f4bdf6c4SBarry Smith     for (i= 0; i< nvars; i++) {
7118b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,yy+(size*i)));
712f4bdf6c4SBarry Smith     }
713f4bdf6c4SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
714f4bdf6c4SBarry Smith   } else {      /* nodal ordering must be mapped to variable ordering for sys_pfmg */
715f4bdf6c4SBarry Smith     PetscScalar *z;
7164ddd07fcSJed Brown     PetscInt    j, k;
717f4bdf6c4SBarry Smith 
718785e854fSJed Brown     ierr = PetscMalloc1(nvars*size,&z);CHKERRQ(ierr);
719fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorSetConstantValues,(mx->ss_b,0.0));
720*b026d285SBarry Smith     ierr = VecGetArrayRead(x,&xx);CHKERRQ(ierr);
721f4bdf6c4SBarry Smith 
722f4bdf6c4SBarry Smith     /* transform nodal to hypre's variable ordering for sys_pfmg */
723f4bdf6c4SBarry Smith     for (i= 0; i< size; i++) {
724f4bdf6c4SBarry Smith       k= i*nvars;
7258865f1eaSKarl Rupp       for (j= 0; j< nvars; j++) z[j*size+i]= xx[k+j];
726f4bdf6c4SBarry Smith     }
727f4bdf6c4SBarry Smith     for (i= 0; i< nvars; i++) {
7288b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructVectorSetBoxValues,(mx->ss_b,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
729f4bdf6c4SBarry Smith     }
730*b026d285SBarry Smith     ierr = VecRestoreArrayRead(x,&xx);CHKERRQ(ierr);
731fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructVectorAssemble,(mx->ss_b));
732fd3f9acdSBarry Smith     PetscStackCallStandard(HYPRE_SStructMatrixMatvec,(1.0,mx->ss_mat,mx->ss_b,0.0,mx->ss_x));
733f4bdf6c4SBarry Smith 
734f4bdf6c4SBarry Smith     /* copy solution values back to PETSc */
735f4bdf6c4SBarry Smith     ierr = VecGetArray(y,&yy);CHKERRQ(ierr);
736f4bdf6c4SBarry Smith     for (i= 0; i< nvars; i++) {
7378b1f7689SBarry Smith       PetscStackCallStandard(HYPRE_SStructVectorGetBoxValues,(mx->ss_x,part,(HYPRE_Int *)ilower,(HYPRE_Int *)iupper,i,z+(size*i)));
738f4bdf6c4SBarry Smith     }
739f4bdf6c4SBarry Smith     /* transform hypre's variable ordering for sys_pfmg to nodal ordering */
740f4bdf6c4SBarry Smith     for (i= 0; i< size; i++) {
741f4bdf6c4SBarry Smith       k= i*nvars;
7428865f1eaSKarl Rupp       for (j= 0; j< nvars; j++) yy[k+j]= z[j*size+i];
743f4bdf6c4SBarry Smith     }
744f4bdf6c4SBarry Smith     ierr = VecRestoreArray(y,&yy);CHKERRQ(ierr);
745f4bdf6c4SBarry Smith     ierr = PetscFree(z);CHKERRQ(ierr);
746f4bdf6c4SBarry Smith   }
747f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
748f4bdf6c4SBarry Smith }
749f4bdf6c4SBarry Smith 
750f4bdf6c4SBarry Smith PetscErrorCode MatAssemblyEnd_HYPRESStruct(Mat mat,MatAssemblyType mode)
751f4bdf6c4SBarry Smith {
752f4bdf6c4SBarry Smith   Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data;
753f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
754f4bdf6c4SBarry Smith 
755f4bdf6c4SBarry Smith   PetscFunctionBegin;
756fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructMatrixAssemble,(ex->ss_mat));
757f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
758f4bdf6c4SBarry Smith }
759f4bdf6c4SBarry Smith 
760f4bdf6c4SBarry Smith PetscErrorCode MatZeroEntries_HYPRESStruct(Mat mat)
761f4bdf6c4SBarry Smith {
762f4bdf6c4SBarry Smith   PetscFunctionBegin;
763f4bdf6c4SBarry Smith   /* before the DMDA is set to the matrix the zero doesn't need to do anything */
764f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
765f4bdf6c4SBarry Smith }
766f4bdf6c4SBarry Smith 
767f4bdf6c4SBarry Smith PetscErrorCode MatDestroy_HYPRESStruct(Mat mat)
768f4bdf6c4SBarry Smith {
769f4bdf6c4SBarry Smith   Mat_HYPRESStruct *ex = (Mat_HYPRESStruct*) mat->data;
770f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
771f4bdf6c4SBarry Smith 
772f4bdf6c4SBarry Smith   PetscFunctionBegin;
773fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructGraphDestroy,(ex->ss_graph));
774fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructMatrixDestroy,(ex->ss_mat));
775fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorDestroy,(ex->ss_x));
776fd3f9acdSBarry Smith   PetscStackCallStandard(HYPRE_SStructVectorDestroy,(ex->ss_b));
777f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
778f4bdf6c4SBarry Smith }
779f4bdf6c4SBarry Smith 
7808cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_HYPRESStruct(Mat B)
781f4bdf6c4SBarry Smith {
782f4bdf6c4SBarry Smith   Mat_HYPRESStruct *ex;
783f4bdf6c4SBarry Smith   PetscErrorCode   ierr;
784f4bdf6c4SBarry Smith 
785f4bdf6c4SBarry Smith   PetscFunctionBegin;
786b00a9115SJed Brown   ierr         = PetscNewLog(B,&ex);CHKERRQ(ierr);
787f4bdf6c4SBarry Smith   B->data      = (void*)ex;
788f4bdf6c4SBarry Smith   B->rmap->bs  = 1;
789f4bdf6c4SBarry Smith   B->assembled = PETSC_FALSE;
790f4bdf6c4SBarry Smith 
791f4bdf6c4SBarry Smith   B->insertmode = NOT_SET_VALUES;
792f4bdf6c4SBarry Smith 
793f4bdf6c4SBarry Smith   B->ops->assemblyend = MatAssemblyEnd_HYPRESStruct;
794f4bdf6c4SBarry Smith   B->ops->mult        = MatMult_HYPRESStruct;
795f4bdf6c4SBarry Smith   B->ops->zeroentries = MatZeroEntries_HYPRESStruct;
796f4bdf6c4SBarry Smith   B->ops->destroy     = MatDestroy_HYPRESStruct;
797f4bdf6c4SBarry Smith 
798f4bdf6c4SBarry Smith   ex->needsinitialization = PETSC_TRUE;
799f4bdf6c4SBarry Smith 
800ce94432eSBarry Smith   ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)B),&(ex->hcomm));CHKERRQ(ierr);
801bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)B,"MatSetupDM_C",MatSetupDM_HYPRESStruct);CHKERRQ(ierr);
802f4bdf6c4SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)B,MATHYPRESSTRUCT);CHKERRQ(ierr);
803f4bdf6c4SBarry Smith   PetscFunctionReturn(0);
804f4bdf6c4SBarry Smith }
805f4bdf6c4SBarry Smith 
806f4bdf6c4SBarry Smith 
807f4bdf6c4SBarry Smith 
808