xref: /petsc/src/dm/impls/da/fdda.c (revision 07475bc16356fc37e8c66fcce1957fb7d8feef24)
147c6ae99SBarry Smith 
2b45d2f2cSJed Brown #include <petsc-private/daimpl.h> /*I      "petscdmda.h"     I*/
3*07475bc1SBarry Smith #include <petscmat.h>
4b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
547c6ae99SBarry Smith 
6e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring *);
7e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring *);
8e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring *);
9e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM,ISColoringType,ISColoring *);
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith /*
1247c6ae99SBarry Smith    For ghost i that may be negative or greater than the upper bound this
1347c6ae99SBarry Smith   maps it into the 0:m-1 range using periodicity
1447c6ae99SBarry Smith */
1547c6ae99SBarry Smith #define SetInRange(i,m) ((i < 0) ? m+i:((i >= m) ? i-m:i))
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith #undef __FUNCT__
18aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills_Private"
19ce308e1dSBarry Smith static PetscErrorCode DMDASetBlockFills_Private(const PetscInt *dfill,PetscInt w,PetscInt **rfill)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   PetscErrorCode ierr;
2247c6ae99SBarry Smith   PetscInt       i,j,nz,*fill;
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   PetscFunctionBegin;
2547c6ae99SBarry Smith   if (!dfill) PetscFunctionReturn(0);
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   /* count number nonzeros */
2847c6ae99SBarry Smith   nz = 0;
2947c6ae99SBarry Smith   for (i=0; i<w; i++) {
3047c6ae99SBarry Smith     for (j=0; j<w; j++) {
3147c6ae99SBarry Smith       if (dfill[w*i+j]) nz++;
3247c6ae99SBarry Smith     }
3347c6ae99SBarry Smith   }
3447c6ae99SBarry Smith   ierr = PetscMalloc((nz + w + 1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
3547c6ae99SBarry Smith   /* construct modified CSR storage of nonzero structure */
36ce308e1dSBarry Smith   /*  fill[0 -- w] marks starts of each row of column indices (and end of last row)
37ce308e1dSBarry Smith    so fill[1] - fill[0] gives number of nonzeros in first row etc */
3847c6ae99SBarry Smith   nz = w + 1;
3947c6ae99SBarry Smith   for (i=0; i<w; i++) {
4047c6ae99SBarry Smith     fill[i] = nz;
4147c6ae99SBarry Smith     for (j=0; j<w; j++) {
4247c6ae99SBarry Smith       if (dfill[w*i+j]) {
4347c6ae99SBarry Smith 	fill[nz] = j;
4447c6ae99SBarry Smith 	nz++;
4547c6ae99SBarry Smith       }
4647c6ae99SBarry Smith     }
4747c6ae99SBarry Smith   }
4847c6ae99SBarry Smith   fill[w] = nz;
4947c6ae99SBarry Smith 
5047c6ae99SBarry Smith   *rfill = fill;
5147c6ae99SBarry Smith   PetscFunctionReturn(0);
5247c6ae99SBarry Smith }
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith #undef __FUNCT__
55aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills"
5647c6ae99SBarry Smith /*@
57aa219208SBarry Smith     DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem
58950540a4SJed Brown     of the matrix returned by DMCreateMatrix().
5947c6ae99SBarry Smith 
60aa219208SBarry Smith     Logically Collective on DMDA
6147c6ae99SBarry Smith 
6247c6ae99SBarry Smith     Input Parameter:
6347c6ae99SBarry Smith +   da - the distributed array
6447c6ae99SBarry Smith .   dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block)
6547c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
6647c6ae99SBarry Smith 
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith     Level: developer
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     Notes: This only makes sense when you are doing multicomponent problems but using the
7147c6ae99SBarry Smith        MPIAIJ matrix format
7247c6ae99SBarry Smith 
7347c6ae99SBarry Smith            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
7447c6ae99SBarry Smith        representing coupling and 0 entries for missing coupling. For example
7547c6ae99SBarry Smith $             dfill[9] = {1, 0, 0,
7647c6ae99SBarry Smith $                         1, 1, 0,
7747c6ae99SBarry Smith $                         0, 1, 1}
7847c6ae99SBarry Smith        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
7947c6ae99SBarry Smith        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
8047c6ae99SBarry Smith        diagonal block).
8147c6ae99SBarry Smith 
82aa219208SBarry Smith      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
8347c6ae99SBarry Smith      can be represented in the dfill, ofill format
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith    Contributed by Glenn Hammond
8647c6ae99SBarry Smith 
878ddb5d8bSBarry Smith .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
8847c6ae99SBarry Smith 
8947c6ae99SBarry Smith @*/
90ce308e1dSBarry Smith PetscErrorCode  DMDASetBlockFills(DM da,const PetscInt *dfill,const PetscInt *ofill)
9147c6ae99SBarry Smith {
9247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
9347c6ae99SBarry Smith   PetscErrorCode ierr;
94ae4f298aSBarry Smith   PetscInt       i,k,cnt = 1;
9547c6ae99SBarry Smith 
9647c6ae99SBarry Smith   PetscFunctionBegin;
97aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr);
98aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr);
99ae4f298aSBarry Smith 
100ae4f298aSBarry Smith   /* ofillcount tracks the columns of ofill that have any nonzero in thems; the value in each location is the number of
101ae4f298aSBarry Smith    columns to the left with any nonzeros in them plus 1 */
102ae4f298aSBarry Smith   ierr = PetscMalloc(dd->w*sizeof(PetscBool),&dd->ofillcols);CHKERRQ(ierr);
103ae4f298aSBarry Smith   ierr = PetscMemzero(dd->ofillcols,dd->w*sizeof(PetscBool));CHKERRQ(ierr);
104ae4f298aSBarry Smith   for (i=0; i<dd->w; i++) {
105ae4f298aSBarry Smith     for (k=dd->ofill[i]; k<dd->ofill[i+1]; k++) dd->ofillcols[dd->ofill[k]] = 1;
106ae4f298aSBarry Smith   }
107ae4f298aSBarry Smith   for (i=0; i<dd->w; i++) {
108ae4f298aSBarry Smith     if (dd->ofillcols[i]) {
109ae4f298aSBarry Smith       dd->ofillcols[i] = cnt++;
110ae4f298aSBarry Smith     }
111ae4f298aSBarry Smith   }
11247c6ae99SBarry Smith   PetscFunctionReturn(0);
11347c6ae99SBarry Smith }
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith 
11647c6ae99SBarry Smith #undef __FUNCT__
117e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA"
11819fd82e9SBarry Smith PetscErrorCode  DMCreateColoring_DA(DM da,ISColoringType ctype,MatType mtype,ISColoring *coloring)
11947c6ae99SBarry Smith {
12047c6ae99SBarry Smith   PetscErrorCode   ierr;
12147c6ae99SBarry Smith   PetscInt         dim,m,n,p,nc;
1221321219cSEthan Coon   DMDABoundaryType bx,by,bz;
12347c6ae99SBarry Smith   MPI_Comm         comm;
12447c6ae99SBarry Smith   PetscMPIInt      size;
12547c6ae99SBarry Smith   PetscBool        isBAIJ;
12647c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith   PetscFunctionBegin;
12947c6ae99SBarry Smith   /*
13047c6ae99SBarry Smith                                   m
13147c6ae99SBarry Smith           ------------------------------------------------------
13247c6ae99SBarry Smith          |                                                     |
13347c6ae99SBarry Smith          |                                                     |
13447c6ae99SBarry Smith          |               ----------------------                |
13547c6ae99SBarry Smith          |               |                    |                |
13647c6ae99SBarry Smith       n  |           yn  |                    |                |
13747c6ae99SBarry Smith          |               |                    |                |
13847c6ae99SBarry Smith          |               .---------------------                |
13947c6ae99SBarry Smith          |             (xs,ys)     xn                          |
14047c6ae99SBarry Smith          |            .                                        |
14147c6ae99SBarry Smith          |         (gxs,gys)                                   |
14247c6ae99SBarry Smith          |                                                     |
14347c6ae99SBarry Smith           -----------------------------------------------------
14447c6ae99SBarry Smith   */
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith   /*
14747c6ae99SBarry Smith          nc - number of components per grid point
14847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
14947c6ae99SBarry Smith 
15047c6ae99SBarry Smith   */
1511321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&bx,&by,&bz,0);CHKERRQ(ierr);
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
15447c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
15547c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED){
15647c6ae99SBarry Smith     if (size == 1) {
15747c6ae99SBarry Smith       ctype = IS_COLORING_GLOBAL;
15847c6ae99SBarry Smith     } else if (dim > 1){
1591321219cSEthan Coon       if ((m==1 && bx == DMDA_BOUNDARY_PERIODIC) || (n==1 && by == DMDA_BOUNDARY_PERIODIC) || (p==1 && bz == DMDA_BOUNDARY_PERIODIC)){
16047c6ae99SBarry Smith         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"IS_COLORING_GHOSTED cannot be used for periodic boundary condition having both ends of the domain  on the same process");
16147c6ae99SBarry Smith       }
16247c6ae99SBarry Smith     }
16347c6ae99SBarry Smith   }
16447c6ae99SBarry Smith 
165aa219208SBarry Smith   /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ
16647c6ae99SBarry Smith      matrices is for the blocks, not the individual matrix elements  */
1674833614aSSatish Balay   ierr = PetscStrcmp(mtype,MATBAIJ,&isBAIJ);CHKERRQ(ierr);
1684833614aSSatish Balay   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);}
16947c6ae99SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);}
17047c6ae99SBarry Smith   if (isBAIJ) {
17147c6ae99SBarry Smith     dd->w = 1;
17247c6ae99SBarry Smith     dd->xs = dd->xs/nc;
17347c6ae99SBarry Smith     dd->xe = dd->xe/nc;
17447c6ae99SBarry Smith     dd->Xs = dd->Xs/nc;
17547c6ae99SBarry Smith     dd->Xe = dd->Xe/nc;
17647c6ae99SBarry Smith   }
17747c6ae99SBarry Smith 
17847c6ae99SBarry Smith   /*
179aa219208SBarry Smith      We do not provide a getcoloring function in the DMDA operations because
180aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
18147c6ae99SBarry Smith    more low-level then matrices.
18247c6ae99SBarry Smith   */
18347c6ae99SBarry Smith   if (dim == 1) {
184e727c939SJed Brown     ierr = DMCreateColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
18547c6ae99SBarry Smith   } else if (dim == 2) {
186e727c939SJed Brown     ierr =  DMCreateColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
18747c6ae99SBarry Smith   } else if (dim == 3) {
188e727c939SJed Brown     ierr =  DMCreateColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
18971cd77b2SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim);
19047c6ae99SBarry Smith   if (isBAIJ) {
19147c6ae99SBarry Smith     dd->w = nc;
19247c6ae99SBarry Smith     dd->xs = dd->xs*nc;
19347c6ae99SBarry Smith     dd->xe = dd->xe*nc;
19447c6ae99SBarry Smith     dd->Xs = dd->Xs*nc;
19547c6ae99SBarry Smith     dd->Xe = dd->Xe*nc;
19647c6ae99SBarry Smith   }
19747c6ae99SBarry Smith   PetscFunctionReturn(0);
19847c6ae99SBarry Smith }
19947c6ae99SBarry Smith 
20047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith #undef __FUNCT__
203e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_MPIAIJ"
204e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
20547c6ae99SBarry Smith {
20647c6ae99SBarry Smith   PetscErrorCode   ierr;
20747c6ae99SBarry Smith   PetscInt         xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col;
20847c6ae99SBarry Smith   PetscInt         ncolors;
20947c6ae99SBarry Smith   MPI_Comm         comm;
2101321219cSEthan Coon   DMDABoundaryType bx,by;
211aa219208SBarry Smith   DMDAStencilType  st;
21247c6ae99SBarry Smith   ISColoringValue  *colors;
21347c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith   PetscFunctionBegin;
21647c6ae99SBarry Smith   /*
21747c6ae99SBarry Smith          nc - number of components per grid point
21847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
21947c6ae99SBarry Smith 
22047c6ae99SBarry Smith   */
2211321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
22247c6ae99SBarry Smith   col    = 2*s + 1;
223aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
224aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
22547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
22647c6ae99SBarry Smith 
22747c6ae99SBarry Smith   /* special case as taught to us by Paul Hovland */
228aa219208SBarry Smith   if (st == DMDA_STENCIL_STAR && s == 1) {
229e727c939SJed Brown     ierr = DMCreateColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
23047c6ae99SBarry Smith   } else {
23147c6ae99SBarry Smith 
23266a15934SBarry Smith     if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\
23347c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", m, col);
23466a15934SBarry Smith     if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\
23547c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", n, col);
23647c6ae99SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
23747c6ae99SBarry Smith       if (!dd->localcoloring) {
23847c6ae99SBarry Smith 	ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
23947c6ae99SBarry Smith 	ii = 0;
24047c6ae99SBarry Smith 	for (j=ys; j<ys+ny; j++) {
24147c6ae99SBarry Smith 	  for (i=xs; i<xs+nx; i++) {
24247c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
24347c6ae99SBarry Smith 	      colors[ii++] = k + nc*((i % col) + col*(j % col));
24447c6ae99SBarry Smith 	    }
24547c6ae99SBarry Smith 	  }
24647c6ae99SBarry Smith 	}
24747c6ae99SBarry Smith         ncolors = nc + nc*(col-1 + col*(col-1));
24847c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
24947c6ae99SBarry Smith       }
25047c6ae99SBarry Smith       *coloring = dd->localcoloring;
25147c6ae99SBarry Smith     } else if (ctype == IS_COLORING_GHOSTED) {
25247c6ae99SBarry Smith       if (!dd->ghostedcoloring) {
25347c6ae99SBarry Smith 	ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
25447c6ae99SBarry Smith 	ii = 0;
25547c6ae99SBarry Smith 	for (j=gys; j<gys+gny; j++) {
25647c6ae99SBarry Smith 	  for (i=gxs; i<gxs+gnx; i++) {
25747c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
25847c6ae99SBarry Smith 	      /* the complicated stuff is to handle periodic boundaries */
25947c6ae99SBarry Smith 	      colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
26047c6ae99SBarry Smith 	    }
26147c6ae99SBarry Smith 	  }
26247c6ae99SBarry Smith 	}
26347c6ae99SBarry Smith         ncolors = nc + nc*(col - 1 + col*(col-1));
26447c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
26547c6ae99SBarry Smith         /* PetscIntView(ncolors,(PetscInt *)colors,0); */
26647c6ae99SBarry Smith 
26747c6ae99SBarry Smith 	ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
26847c6ae99SBarry Smith       }
26947c6ae99SBarry Smith       *coloring = dd->ghostedcoloring;
27047c6ae99SBarry Smith     } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
27147c6ae99SBarry Smith   }
27247c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
27347c6ae99SBarry Smith   PetscFunctionReturn(0);
27447c6ae99SBarry Smith }
27547c6ae99SBarry Smith 
27647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
27747c6ae99SBarry Smith 
27847c6ae99SBarry Smith #undef __FUNCT__
279e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_3d_MPIAIJ"
280e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
28147c6ae99SBarry Smith {
28247c6ae99SBarry Smith   PetscErrorCode    ierr;
28347c6ae99SBarry Smith   PetscInt          xs,ys,nx,ny,i,j,gxs,gys,gnx,gny,m,n,p,dim,s,k,nc,col,zs,gzs,ii,l,nz,gnz,M,N,P;
28447c6ae99SBarry Smith   PetscInt          ncolors;
28547c6ae99SBarry Smith   MPI_Comm          comm;
2861321219cSEthan Coon   DMDABoundaryType  bx,by,bz;
287aa219208SBarry Smith   DMDAStencilType   st;
28847c6ae99SBarry Smith   ISColoringValue   *colors;
28947c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
29047c6ae99SBarry Smith 
29147c6ae99SBarry Smith   PetscFunctionBegin;
29247c6ae99SBarry Smith   /*
29347c6ae99SBarry Smith          nc - number of components per grid point
29447c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
29547c6ae99SBarry Smith 
29647c6ae99SBarry Smith   */
2971321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
29847c6ae99SBarry Smith   col    = 2*s + 1;
29966a15934SBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
30047c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
30166a15934SBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
30247c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
30366a15934SBarry Smith   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
30447c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
30547c6ae99SBarry Smith 
306aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
307aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
30847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
30947c6ae99SBarry Smith 
31047c6ae99SBarry Smith   /* create the coloring */
31147c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
31247c6ae99SBarry Smith     if (!dd->localcoloring) {
31347c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
31447c6ae99SBarry Smith       ii = 0;
31547c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
31647c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
31747c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
31847c6ae99SBarry Smith             for (l=0; l<nc; l++) {
31947c6ae99SBarry Smith               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
32047c6ae99SBarry Smith             }
32147c6ae99SBarry Smith           }
32247c6ae99SBarry Smith         }
32347c6ae99SBarry Smith       }
32447c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
32547c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr);
32647c6ae99SBarry Smith     }
32747c6ae99SBarry Smith     *coloring = dd->localcoloring;
32847c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
32947c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
33047c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
33147c6ae99SBarry Smith       ii = 0;
33247c6ae99SBarry Smith       for (k=gzs; k<gzs+gnz; k++) {
33347c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
33447c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
33547c6ae99SBarry Smith             for (l=0; l<nc; l++) {
33647c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
33747c6ae99SBarry Smith               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
33847c6ae99SBarry Smith             }
33947c6ae99SBarry Smith           }
34047c6ae99SBarry Smith         }
34147c6ae99SBarry Smith       }
34247c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
34347c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
34447c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
34547c6ae99SBarry Smith     }
34647c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
34747c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
34847c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
34947c6ae99SBarry Smith   PetscFunctionReturn(0);
35047c6ae99SBarry Smith }
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
35347c6ae99SBarry Smith 
35447c6ae99SBarry Smith #undef __FUNCT__
355e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_1d_MPIAIJ"
356e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
35747c6ae99SBarry Smith {
35847c6ae99SBarry Smith   PetscErrorCode    ierr;
35947c6ae99SBarry Smith   PetscInt          xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
36047c6ae99SBarry Smith   PetscInt          ncolors;
36147c6ae99SBarry Smith   MPI_Comm          comm;
3621321219cSEthan Coon   DMDABoundaryType  bx;
36347c6ae99SBarry Smith   ISColoringValue   *colors;
36447c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
36547c6ae99SBarry Smith 
36647c6ae99SBarry Smith   PetscFunctionBegin;
36747c6ae99SBarry Smith   /*
36847c6ae99SBarry Smith          nc - number of components per grid point
36947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
37047c6ae99SBarry Smith 
37147c6ae99SBarry Smith   */
3721321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
37347c6ae99SBarry Smith   col    = 2*s + 1;
37447c6ae99SBarry Smith 
37566a15934SBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points %d is divisible\n\
37631e6f798SBarry Smith                                                           by 2*stencil_width + 1 %d\n",(int)m,(int)col);
37747c6ae99SBarry Smith 
378aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
379aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
38047c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
38147c6ae99SBarry Smith 
38247c6ae99SBarry Smith   /* create the coloring */
38347c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
38447c6ae99SBarry Smith     if (!dd->localcoloring) {
38547c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
386ae4f298aSBarry Smith       if (dd->ofillcols) {
387ae4f298aSBarry Smith         PetscInt tc = 0;
388ae4f298aSBarry Smith         for (i=0; i<nc; i++) tc += (PetscInt) (dd->ofillcols[i] > 0);
389ae4f298aSBarry Smith         i1 = 0;
390ae4f298aSBarry Smith         for (i=xs; i<xs+nx; i++) {
391ae4f298aSBarry Smith           for (l=0; l<nc; l++) {
392ae4f298aSBarry Smith             if (dd->ofillcols[l] && (i % col)) {
393ae4f298aSBarry Smith               colors[i1++] =  nc - 1 + tc*((i % col) - 1) + dd->ofillcols[l];
394ae4f298aSBarry Smith             } else {
395ae4f298aSBarry Smith               colors[i1++] = l;
396ae4f298aSBarry Smith             }
397ae4f298aSBarry Smith           }
398ae4f298aSBarry Smith         }
399ae4f298aSBarry Smith         ncolors = nc + 2*s*tc;
400ae4f298aSBarry Smith       } else {
40147c6ae99SBarry Smith         i1 = 0;
40247c6ae99SBarry Smith         for (i=xs; i<xs+nx; i++) {
40347c6ae99SBarry Smith           for (l=0; l<nc; l++) {
40447c6ae99SBarry Smith             colors[i1++] = l + nc*(i % col);
40547c6ae99SBarry Smith           }
40647c6ae99SBarry Smith         }
40747c6ae99SBarry Smith         ncolors = nc + nc*(col-1);
408ae4f298aSBarry Smith       }
40947c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr);
41047c6ae99SBarry Smith     }
41147c6ae99SBarry Smith     *coloring = dd->localcoloring;
41247c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
41347c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
41447c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
41547c6ae99SBarry Smith       i1 = 0;
41647c6ae99SBarry Smith       for (i=gxs; i<gxs+gnx; i++) {
41747c6ae99SBarry Smith         for (l=0; l<nc; l++) {
41847c6ae99SBarry Smith           /* the complicated stuff is to handle periodic boundaries */
41947c6ae99SBarry Smith           colors[i1++] = l + nc*(SetInRange(i,m) % col);
42047c6ae99SBarry Smith         }
42147c6ae99SBarry Smith       }
42247c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
42347c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
42447c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
42547c6ae99SBarry Smith     }
42647c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
42747c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
42847c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
42947c6ae99SBarry Smith   PetscFunctionReturn(0);
43047c6ae99SBarry Smith }
43147c6ae99SBarry Smith 
43247c6ae99SBarry Smith #undef __FUNCT__
433e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_5pt_MPIAIJ"
434e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
43547c6ae99SBarry Smith {
43647c6ae99SBarry Smith   PetscErrorCode    ierr;
43747c6ae99SBarry Smith   PetscInt          xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
43847c6ae99SBarry Smith   PetscInt          ncolors;
43947c6ae99SBarry Smith   MPI_Comm          comm;
4401321219cSEthan Coon   DMDABoundaryType  bx,by;
44147c6ae99SBarry Smith   ISColoringValue   *colors;
44247c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
44347c6ae99SBarry Smith 
44447c6ae99SBarry Smith   PetscFunctionBegin;
44547c6ae99SBarry Smith   /*
44647c6ae99SBarry Smith          nc - number of components per grid point
44747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
44847c6ae99SBarry Smith 
44947c6ae99SBarry Smith   */
4501321219cSEthan Coon   ierr   = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);CHKERRQ(ierr);
451aa219208SBarry Smith   ierr   = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
452aa219208SBarry Smith   ierr   = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
45347c6ae99SBarry Smith   ierr   = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
45447c6ae99SBarry Smith 
45566a15934SBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % 5)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible by 5\n");
45666a15934SBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % 5)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible by 5\n");
45747c6ae99SBarry Smith 
45847c6ae99SBarry Smith   /* create the coloring */
45947c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
46047c6ae99SBarry Smith     if (!dd->localcoloring) {
46147c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
46247c6ae99SBarry Smith       ii = 0;
46347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
46447c6ae99SBarry Smith 	for (i=xs; i<xs+nx; i++) {
46547c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
46647c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*j+i) % 5);
46747c6ae99SBarry Smith 	  }
46847c6ae99SBarry Smith 	}
46947c6ae99SBarry Smith       }
47047c6ae99SBarry Smith       ncolors = 5*nc;
47147c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
47247c6ae99SBarry Smith     }
47347c6ae99SBarry Smith     *coloring = dd->localcoloring;
47447c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
47547c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
47647c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
47747c6ae99SBarry Smith       ii = 0;
47847c6ae99SBarry Smith       for (j=gys; j<gys+gny; j++) {
47947c6ae99SBarry Smith 	for (i=gxs; i<gxs+gnx; i++) {
48047c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
48147c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
48247c6ae99SBarry Smith 	  }
48347c6ae99SBarry Smith 	}
48447c6ae99SBarry Smith       }
48547c6ae99SBarry Smith       ncolors = 5*nc;
48647c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
48747c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
48847c6ae99SBarry Smith     }
48947c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
49047c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
49147c6ae99SBarry Smith   PetscFunctionReturn(0);
49247c6ae99SBarry Smith }
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith /* =========================================================================== */
495950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM,Mat);
496ce308e1dSBarry Smith extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM,Mat);
497950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM,Mat);
498950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM,Mat);
499950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM,Mat);
500950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM,Mat);
501950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM,Mat);
502950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM,Mat);
503950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM,Mat);
504950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM,Mat);
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith #undef __FUNCT__
507c688c046SMatthew G Knepley #define __FUNCT__ "MatSetupDM"
5088bbdbebaSMatthew G Knepley /*@C
509c688c046SMatthew G Knepley    MatSetupDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith    Logically Collective on Mat
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith    Input Parameters:
51447c6ae99SBarry Smith +  mat - the matrix
51547c6ae99SBarry Smith -  da - the da
51647c6ae99SBarry Smith 
51747c6ae99SBarry Smith    Level: intermediate
51847c6ae99SBarry Smith 
51947c6ae99SBarry Smith @*/
520c688c046SMatthew G Knepley PetscErrorCode MatSetupDM(Mat mat,DM da)
52147c6ae99SBarry Smith {
52247c6ae99SBarry Smith   PetscErrorCode ierr;
52347c6ae99SBarry Smith 
52447c6ae99SBarry Smith   PetscFunctionBegin;
52547c6ae99SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
52647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
527c688c046SMatthew G Knepley   ierr = PetscTryMethod(mat,"MatSetupDM_C",(Mat,DM),(mat,da));CHKERRQ(ierr);
52847c6ae99SBarry Smith   PetscFunctionReturn(0);
52947c6ae99SBarry Smith }
53047c6ae99SBarry Smith 
53147c6ae99SBarry Smith EXTERN_C_BEGIN
53247c6ae99SBarry Smith #undef __FUNCT__
53347c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA"
5347087cfbeSBarry Smith PetscErrorCode  MatView_MPI_DA(Mat A,PetscViewer viewer)
53547c6ae99SBarry Smith {
5369a42bb27SBarry Smith   DM             da;
53747c6ae99SBarry Smith   PetscErrorCode ierr;
53847c6ae99SBarry Smith   const char     *prefix;
53947c6ae99SBarry Smith   Mat            Anatural;
54047c6ae99SBarry Smith   AO             ao;
54147c6ae99SBarry Smith   PetscInt       rstart,rend,*petsc,i;
54247c6ae99SBarry Smith   IS             is;
54347c6ae99SBarry Smith   MPI_Comm       comm;
54474388724SJed Brown   PetscViewerFormat format;
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith   PetscFunctionBegin;
54774388724SJed Brown   /* Check whether we are just printing info, in which case MatView() already viewed everything we wanted to view */
54874388724SJed Brown   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
54974388724SJed Brown   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
55074388724SJed Brown 
55147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
552c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
553aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
55447c6ae99SBarry Smith 
555aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
55647c6ae99SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
55747c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr);
55847c6ae99SBarry Smith   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
55947c6ae99SBarry Smith   ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr);
56047c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith   /* call viewer on natural ordering */
56347c6ae99SBarry Smith   ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr);
564fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
56547c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr);
56647c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr);
56747c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr);
56847c6ae99SBarry Smith   ierr = MatView(Anatural,viewer);CHKERRQ(ierr);
569fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
57047c6ae99SBarry Smith   PetscFunctionReturn(0);
57147c6ae99SBarry Smith }
57247c6ae99SBarry Smith EXTERN_C_END
57347c6ae99SBarry Smith 
57447c6ae99SBarry Smith EXTERN_C_BEGIN
57547c6ae99SBarry Smith #undef __FUNCT__
57647c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA"
5777087cfbeSBarry Smith PetscErrorCode  MatLoad_MPI_DA(Mat A,PetscViewer viewer)
57847c6ae99SBarry Smith {
5799a42bb27SBarry Smith   DM             da;
58047c6ae99SBarry Smith   PetscErrorCode ierr;
58147c6ae99SBarry Smith   Mat            Anatural,Aapp;
58247c6ae99SBarry Smith   AO             ao;
58347c6ae99SBarry Smith   PetscInt       rstart,rend,*app,i;
58447c6ae99SBarry Smith   IS             is;
58547c6ae99SBarry Smith   MPI_Comm       comm;
58647c6ae99SBarry Smith 
58747c6ae99SBarry Smith   PetscFunctionBegin;
58847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
589c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
590aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
59147c6ae99SBarry Smith 
59247c6ae99SBarry Smith   /* Load the matrix in natural ordering */
59347c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,&Anatural);CHKERRQ(ierr);
59447c6ae99SBarry Smith   ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr);
59547c6ae99SBarry Smith   ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
59647c6ae99SBarry Smith   ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr);
59747c6ae99SBarry Smith 
59847c6ae99SBarry Smith   /* Map natural ordering to application ordering and create IS */
599aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
60047c6ae99SBarry Smith   ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr);
60147c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr);
60247c6ae99SBarry Smith   for (i=rstart; i<rend; i++) app[i-rstart] = i;
60347c6ae99SBarry Smith   ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr);
60447c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
60547c6ae99SBarry Smith 
60647c6ae99SBarry Smith   /* Do permutation and replace header */
60747c6ae99SBarry Smith   ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr);
60847c6ae99SBarry Smith   ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr);
609fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
610fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
61147c6ae99SBarry Smith   PetscFunctionReturn(0);
61247c6ae99SBarry Smith }
61347c6ae99SBarry Smith EXTERN_C_END
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith #undef __FUNCT__
616950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA"
61719fd82e9SBarry Smith PetscErrorCode DMCreateMatrix_DA(DM da, MatType mtype,Mat *J)
61847c6ae99SBarry Smith {
61947c6ae99SBarry Smith   PetscErrorCode ierr;
62047c6ae99SBarry Smith   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
62147c6ae99SBarry Smith   Mat            A;
62247c6ae99SBarry Smith   MPI_Comm       comm;
62319fd82e9SBarry Smith   MatType        Atype;
62437d0c07bSMatthew G Knepley   PetscSection   section, sectionGlobal;
62547c6ae99SBarry Smith   void           (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL;
62647c6ae99SBarry Smith   MatType        ttype[256];
62747c6ae99SBarry Smith   PetscBool      flg;
62847c6ae99SBarry Smith   PetscMPIInt    size;
62947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
63047c6ae99SBarry Smith 
63147c6ae99SBarry Smith   PetscFunctionBegin;
63247c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
63347c6ae99SBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
63447c6ae99SBarry Smith #endif
6355da5aae0SJed Brown   if (!mtype) mtype = MATAIJ;
63647c6ae99SBarry Smith   ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr);
637aa219208SBarry Smith   ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DMDA options","Mat");CHKERRQ(ierr);
638dd85299cSBarry Smith   ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr);
63947c6ae99SBarry Smith   ierr = PetscOptionsEnd();
64047c6ae99SBarry Smith 
64137d0c07bSMatthew G Knepley   ierr = DMGetDefaultSection(da, &section);CHKERRQ(ierr);
64237d0c07bSMatthew G Knepley   if (section) {
64337d0c07bSMatthew G Knepley     PetscInt  bs = -1;
64437d0c07bSMatthew G Knepley     PetscInt  localSize;
64537d0c07bSMatthew G Knepley     PetscBool isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isSymmetric;
64637d0c07bSMatthew G Knepley 
64737d0c07bSMatthew G Knepley     ierr = DMGetDefaultGlobalSection(da, &sectionGlobal);CHKERRQ(ierr);
64837d0c07bSMatthew G Knepley     ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize);CHKERRQ(ierr);
64937d0c07bSMatthew G Knepley     ierr = MatCreate(((PetscObject) da)->comm, J);CHKERRQ(ierr);
65037d0c07bSMatthew G Knepley     ierr = MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
65137d0c07bSMatthew G Knepley     ierr = MatSetType(*J, mtype);CHKERRQ(ierr);
65237d0c07bSMatthew G Knepley     ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
65337d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSHELL, &isShell);CHKERRQ(ierr);
65437d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATBAIJ, &isBlock);CHKERRQ(ierr);
65537d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSEQBAIJ, &isSeqBlock);CHKERRQ(ierr);
65637d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATMPIBAIJ, &isMPIBlock);CHKERRQ(ierr);
65737d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr);
65837d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr);
65937d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr);
66037d0c07bSMatthew G Knepley     /* Check for symmetric storage */
66137d0c07bSMatthew G Knepley     isSymmetric = (PetscBool) (isSymBlock || isSymSeqBlock || isSymMPIBlock);
66237d0c07bSMatthew G Knepley     if (isSymmetric) {
66337d0c07bSMatthew G Knepley       ierr = MatSetOption(*J, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);
66437d0c07bSMatthew G Knepley     }
66537d0c07bSMatthew G Knepley     if (!isShell) {
6664020307fSSatish Balay       /* PetscBool fillMatrix = (PetscBool) !da->prealloc_only; */
66737d0c07bSMatthew G Knepley       PetscInt *dnz, *onz, *dnzu, *onzu, bsLocal;
66837d0c07bSMatthew G Knepley 
66937d0c07bSMatthew G Knepley       if (bs < 0) {
67037d0c07bSMatthew G Knepley         if (isBlock || isSeqBlock || isMPIBlock || isSymBlock || isSymSeqBlock || isSymMPIBlock) {
67137d0c07bSMatthew G Knepley           PetscInt pStart, pEnd, p, dof;
67237d0c07bSMatthew G Knepley 
67337d0c07bSMatthew G Knepley           ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
67437d0c07bSMatthew G Knepley           for (p = pStart; p < pEnd; ++p) {
67537d0c07bSMatthew G Knepley             ierr = PetscSectionGetDof(sectionGlobal, p, &dof);CHKERRQ(ierr);
67637d0c07bSMatthew G Knepley             if (dof) {
67737d0c07bSMatthew G Knepley               bs = dof;
67837d0c07bSMatthew G Knepley               break;
67937d0c07bSMatthew G Knepley             }
68037d0c07bSMatthew G Knepley           }
68137d0c07bSMatthew G Knepley         } else {
68237d0c07bSMatthew G Knepley           bs = 1;
68337d0c07bSMatthew G Knepley         }
68437d0c07bSMatthew G Knepley         /* Must have same blocksize on all procs (some might have no points) */
68537d0c07bSMatthew G Knepley         bsLocal = bs;
68637d0c07bSMatthew G Knepley         ierr = MPI_Allreduce(&bsLocal, &bs, 1, MPIU_INT, MPI_MAX, ((PetscObject) da)->comm);CHKERRQ(ierr);
68737d0c07bSMatthew G Knepley       }
68837d0c07bSMatthew G Knepley       ierr = PetscMalloc4(localSize/bs, PetscInt, &dnz, localSize/bs, PetscInt, &onz, localSize/bs, PetscInt, &dnzu, localSize/bs, PetscInt, &onzu);CHKERRQ(ierr);
68937d0c07bSMatthew G Knepley       ierr = PetscMemzero(dnz,  localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
69037d0c07bSMatthew G Knepley       ierr = PetscMemzero(onz,  localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
69137d0c07bSMatthew G Knepley       ierr = PetscMemzero(dnzu, localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
69237d0c07bSMatthew G Knepley       ierr = PetscMemzero(onzu, localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
693552f7358SJed Brown       /* ierr = DMPlexPreallocateOperator(dm, bs, section, sectionGlobal, dnz, onz, dnzu, onzu, *J, fillMatrix);CHKERRQ(ierr); */
69437d0c07bSMatthew G Knepley       ierr = PetscFree4(dnz, onz, dnzu, onzu);CHKERRQ(ierr);
69537d0c07bSMatthew G Knepley     }
69637d0c07bSMatthew G Knepley   }
69747c6ae99SBarry Smith   /*
69847c6ae99SBarry Smith                                   m
69947c6ae99SBarry Smith           ------------------------------------------------------
70047c6ae99SBarry Smith          |                                                     |
70147c6ae99SBarry Smith          |                                                     |
70247c6ae99SBarry Smith          |               ----------------------                |
70347c6ae99SBarry Smith          |               |                    |                |
70447c6ae99SBarry Smith       n  |           ny  |                    |                |
70547c6ae99SBarry Smith          |               |                    |                |
70647c6ae99SBarry Smith          |               .---------------------                |
70747c6ae99SBarry Smith          |             (xs,ys)     nx                          |
70847c6ae99SBarry Smith          |            .                                        |
70947c6ae99SBarry Smith          |         (gxs,gys)                                   |
71047c6ae99SBarry Smith          |                                                     |
71147c6ae99SBarry Smith           -----------------------------------------------------
71247c6ae99SBarry Smith   */
71347c6ae99SBarry Smith 
71447c6ae99SBarry Smith   /*
71547c6ae99SBarry Smith          nc - number of components per grid point
71647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
71747c6ae99SBarry Smith 
71847c6ae99SBarry Smith   */
719e30e807fSPeter Brune   M = dd->M;
720e30e807fSPeter Brune   N = dd->N;
721e30e807fSPeter Brune   P = dd->P;
722e30e807fSPeter Brune   dim = dd->dim;
723e30e807fSPeter Brune   dof = dd->w;
724e30e807fSPeter Brune   /* ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); */
725aa219208SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
72647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
72747c6ae99SBarry Smith   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
72847c6ae99SBarry Smith   ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr);
72919fd82e9SBarry Smith   ierr = MatSetType(A,(MatType)ttype);CHKERRQ(ierr);
73095ee5b0eSBarry Smith   ierr = MatSetDM(A,da);CHKERRQ(ierr);
73147c6ae99SBarry Smith   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
73247c6ae99SBarry Smith   ierr = MatGetType(A,&Atype);CHKERRQ(ierr);
73347c6ae99SBarry Smith   /*
734aa219208SBarry Smith      We do not provide a getmatrix function in the DMDA operations because
735aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
73647c6ae99SBarry Smith    more low-level than matrices. This is kind of cheating but, cause sometimes
737aa219208SBarry Smith    we think of DMDA has higher level than matrices.
73847c6ae99SBarry Smith 
73947c6ae99SBarry Smith      We could switch based on Atype (or mtype), but we do not since the
74047c6ae99SBarry Smith    specialized setting routines depend only the particular preallocation
74147c6ae99SBarry Smith    details of the matrix, not the type itself.
74247c6ae99SBarry Smith   */
74347c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
74447c6ae99SBarry Smith   if (!aij) {
74547c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
74647c6ae99SBarry Smith   }
74747c6ae99SBarry Smith   if (!aij) {
74847c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
74947c6ae99SBarry Smith     if (!baij) {
75047c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
75147c6ae99SBarry Smith     }
75247c6ae99SBarry Smith     if (!baij){
75347c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
75447c6ae99SBarry Smith       if (!sbaij) {
75547c6ae99SBarry Smith         ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
75647c6ae99SBarry Smith       }
75747c6ae99SBarry Smith     }
75847c6ae99SBarry Smith   }
75947c6ae99SBarry Smith   if (aij) {
76047c6ae99SBarry Smith     if (dim == 1) {
761ce308e1dSBarry Smith       if (dd->ofill) {
762ce308e1dSBarry Smith         ierr = DMCreateMatrix_DA_1d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
763ce308e1dSBarry Smith       } else {
764950540a4SJed Brown         ierr = DMCreateMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr);
765ce308e1dSBarry Smith       }
76647c6ae99SBarry Smith     } else if (dim == 2) {
76747c6ae99SBarry Smith       if (dd->ofill) {
768950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
76947c6ae99SBarry Smith       } else {
770950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr);
77147c6ae99SBarry Smith       }
77247c6ae99SBarry Smith     } else if (dim == 3) {
77347c6ae99SBarry Smith       if (dd->ofill) {
774950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
77547c6ae99SBarry Smith       } else {
776950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr);
77747c6ae99SBarry Smith       }
77847c6ae99SBarry Smith     }
77947c6ae99SBarry Smith   } else if (baij) {
78047c6ae99SBarry Smith     if (dim == 2) {
781950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr);
78247c6ae99SBarry Smith     } else if (dim == 3) {
783950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr);
78466a15934SBarry Smith     } else  SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension! Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim);
78547c6ae99SBarry Smith   } else if (sbaij) {
78647c6ae99SBarry Smith     if (dim == 2) {
787950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr);
78847c6ae99SBarry Smith     } else if (dim == 3) {
789950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr);
79066a15934SBarry Smith     } else SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension! Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim);
791869776cdSLisandro Dalcin   } else {
792869776cdSLisandro Dalcin     ISLocalToGlobalMapping ltog,ltogb;
793869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
794869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
7952949035bSJed Brown     ierr = MatSetUp(A);CHKERRQ(ierr);
796869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMapping(A,ltog,ltog);CHKERRQ(ierr);
797869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMappingBlock(A,ltogb,ltogb);CHKERRQ(ierr);
79847c6ae99SBarry Smith   }
799aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
80047c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
801c688c046SMatthew G Knepley   ierr = MatSetDM(A,da);CHKERRQ(ierr);
80247c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
80347c6ae99SBarry Smith   if (size > 1) {
80447c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
80547c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);CHKERRQ(ierr);
80647c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);CHKERRQ(ierr);
80747c6ae99SBarry Smith   }
80847c6ae99SBarry Smith   *J = A;
80947c6ae99SBarry Smith   PetscFunctionReturn(0);
81047c6ae99SBarry Smith }
81147c6ae99SBarry Smith 
81247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
81347c6ae99SBarry Smith #undef __FUNCT__
814950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ"
815950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM da,Mat J)
81647c6ae99SBarry Smith {
81747c6ae99SBarry Smith   PetscErrorCode         ierr;
81847c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny,m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p;
81947c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
82047c6ae99SBarry Smith   MPI_Comm               comm;
82147c6ae99SBarry Smith   PetscScalar            *values;
8221321219cSEthan Coon   DMDABoundaryType       bx,by;
82347c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
824aa219208SBarry Smith   DMDAStencilType        st;
82547c6ae99SBarry Smith 
82647c6ae99SBarry Smith   PetscFunctionBegin;
82747c6ae99SBarry Smith   /*
82847c6ae99SBarry Smith          nc - number of components per grid point
82947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
83047c6ae99SBarry Smith 
83147c6ae99SBarry Smith   */
8321321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
83347c6ae99SBarry Smith   col = 2*s + 1;
834aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
835aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
83647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
83747c6ae99SBarry Smith 
83847c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
8391411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
8401411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
84147c6ae99SBarry Smith 
84247c6ae99SBarry Smith   /* determine the matrix preallocation information */
84347c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
84447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
84547c6ae99SBarry Smith 
8461321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8471321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
84847c6ae99SBarry Smith 
84947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
85047c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
85147c6ae99SBarry Smith 
8521321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8531321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
85447c6ae99SBarry Smith 
85547c6ae99SBarry Smith       cnt  = 0;
85647c6ae99SBarry Smith       for (k=0; k<nc; k++) {
85747c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
85847c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
859aa219208SBarry Smith 	    if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
86047c6ae99SBarry Smith 	      cols[cnt++]  = k + nc*(slot + gnx*l + p);
86147c6ae99SBarry Smith 	    }
86247c6ae99SBarry Smith 	  }
86347c6ae99SBarry Smith 	}
86447c6ae99SBarry Smith 	rows[k] = k + nc*(slot);
86547c6ae99SBarry Smith       }
866784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
86747c6ae99SBarry Smith     }
86847c6ae99SBarry Smith   }
869f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
87047c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
87147c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
87247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
87347c6ae99SBarry Smith 
874784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
875784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
87647c6ae99SBarry Smith 
87747c6ae99SBarry Smith   /*
87847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
87947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
88047c6ae99SBarry Smith     PETSc ordering.
88147c6ae99SBarry Smith   */
882fcfd50ebSBarry Smith   if (!da->prealloc_only) {
88347c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
88447c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
88547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
88647c6ae99SBarry Smith 
8871321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8881321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
88947c6ae99SBarry Smith 
89047c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
89147c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
89247c6ae99SBarry Smith 
8931321219cSEthan Coon 	lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8941321219cSEthan Coon 	lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
89547c6ae99SBarry Smith 
89647c6ae99SBarry Smith 	cnt  = 0;
89747c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
89847c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
89947c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
900aa219208SBarry Smith 	      if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
90147c6ae99SBarry Smith 		cols[cnt++]  = k + nc*(slot + gnx*l + p);
90247c6ae99SBarry Smith 	      }
90347c6ae99SBarry Smith 	    }
90447c6ae99SBarry Smith 	  }
90547c6ae99SBarry Smith 	  rows[k]      = k + nc*(slot);
90647c6ae99SBarry Smith 	}
90747c6ae99SBarry Smith 	ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
90847c6ae99SBarry Smith       }
90947c6ae99SBarry Smith     }
91047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
91147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91347c6ae99SBarry Smith   }
91447c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
91547c6ae99SBarry Smith   PetscFunctionReturn(0);
91647c6ae99SBarry Smith }
91747c6ae99SBarry Smith 
91847c6ae99SBarry Smith #undef __FUNCT__
919950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ_Fill"
920950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J)
92147c6ae99SBarry Smith {
92247c6ae99SBarry Smith   PetscErrorCode         ierr;
92347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
92447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p;
92547c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
92647c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
92747c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
92847c6ae99SBarry Smith   MPI_Comm               comm;
92947c6ae99SBarry Smith   PetscScalar            *values;
9301321219cSEthan Coon   DMDABoundaryType       bx,by;
93147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
932aa219208SBarry Smith   DMDAStencilType        st;
93347c6ae99SBarry Smith 
93447c6ae99SBarry Smith   PetscFunctionBegin;
93547c6ae99SBarry Smith   /*
93647c6ae99SBarry Smith          nc - number of components per grid point
93747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
93847c6ae99SBarry Smith 
93947c6ae99SBarry Smith   */
9401321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
94147c6ae99SBarry Smith   col = 2*s + 1;
942aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
943aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
94447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
94547c6ae99SBarry Smith 
94647c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
9471411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
9481411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
94947c6ae99SBarry Smith 
95047c6ae99SBarry Smith   /* determine the matrix preallocation information */
95147c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
95247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
95347c6ae99SBarry Smith 
9541321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9551321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
95647c6ae99SBarry Smith 
95747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
95847c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
95947c6ae99SBarry Smith 
9601321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9611321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
96247c6ae99SBarry Smith 
96347c6ae99SBarry Smith       for (k=0; k<nc; k++) {
96447c6ae99SBarry Smith         cnt  = 0;
96547c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
96647c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
96747c6ae99SBarry Smith             if (l || p) {
968aa219208SBarry Smith 	      if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
96947c6ae99SBarry Smith                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
97047c6ae99SBarry Smith 		  cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
97147c6ae99SBarry Smith 	      }
97247c6ae99SBarry Smith             } else {
97347c6ae99SBarry Smith 	      if (dfill) {
97447c6ae99SBarry Smith 		for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
97547c6ae99SBarry Smith 		  cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
97647c6ae99SBarry Smith 	      } else {
97747c6ae99SBarry Smith 		for (ifill_col=0; ifill_col<nc; ifill_col++)
97847c6ae99SBarry Smith 		  cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
97947c6ae99SBarry Smith 	      }
98047c6ae99SBarry Smith             }
98147c6ae99SBarry Smith 	  }
98247c6ae99SBarry Smith 	}
98347c6ae99SBarry Smith 	row = k + nc*(slot);
984784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
98547c6ae99SBarry Smith       }
98647c6ae99SBarry Smith     }
98747c6ae99SBarry Smith   }
98847c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
98947c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
99047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
991784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
992784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
99347c6ae99SBarry Smith 
99447c6ae99SBarry Smith   /*
99547c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
99647c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
99747c6ae99SBarry Smith     PETSc ordering.
99847c6ae99SBarry Smith   */
999fcfd50ebSBarry Smith   if (!da->prealloc_only) {
100047c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
100147c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
100247c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
100347c6ae99SBarry Smith 
10041321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10051321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
100647c6ae99SBarry Smith 
100747c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
100847c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
100947c6ae99SBarry Smith 
10101321219cSEthan Coon 	lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10111321219cSEthan Coon 	lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
101247c6ae99SBarry Smith 
101347c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
101447c6ae99SBarry Smith 	  cnt  = 0;
101547c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
101647c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
101747c6ae99SBarry Smith 	      if (l || p) {
1018aa219208SBarry Smith 		if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
101947c6ae99SBarry Smith 		  for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
102047c6ae99SBarry Smith 		    cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
102147c6ae99SBarry Smith 		}
102247c6ae99SBarry Smith 	      } else {
102347c6ae99SBarry Smith 		if (dfill) {
102447c6ae99SBarry Smith 		  for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
102547c6ae99SBarry Smith 		    cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
102647c6ae99SBarry Smith 		} else {
102747c6ae99SBarry Smith 		  for (ifill_col=0; ifill_col<nc; ifill_col++)
102847c6ae99SBarry Smith 		    cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
102947c6ae99SBarry Smith 		}
103047c6ae99SBarry Smith 	      }
103147c6ae99SBarry Smith 	    }
103247c6ae99SBarry Smith 	  }
103347c6ae99SBarry Smith 	  row  = k + nc*(slot);
103447c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
103547c6ae99SBarry Smith 	}
103647c6ae99SBarry Smith       }
103747c6ae99SBarry Smith     }
103847c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
103947c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
104047c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
104147c6ae99SBarry Smith   }
104247c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
104347c6ae99SBarry Smith   PetscFunctionReturn(0);
104447c6ae99SBarry Smith }
104547c6ae99SBarry Smith 
104647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
104747c6ae99SBarry Smith 
104847c6ae99SBarry Smith #undef __FUNCT__
1049950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ"
1050950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM da,Mat J)
105147c6ae99SBarry Smith {
105247c6ae99SBarry Smith   PetscErrorCode         ierr;
105347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
105447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL;
105547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
105647c6ae99SBarry Smith   MPI_Comm               comm;
105747c6ae99SBarry Smith   PetscScalar            *values;
10581321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
105947c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1060aa219208SBarry Smith   DMDAStencilType        st;
106147c6ae99SBarry Smith 
106247c6ae99SBarry Smith   PetscFunctionBegin;
106347c6ae99SBarry Smith   /*
106447c6ae99SBarry Smith          nc - number of components per grid point
106547c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
106647c6ae99SBarry Smith 
106747c6ae99SBarry Smith   */
10681321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
106947c6ae99SBarry Smith   col    = 2*s + 1;
107047c6ae99SBarry Smith 
1071aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1072aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
107347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
107447c6ae99SBarry Smith 
107547c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
10761411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
10771411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
107847c6ae99SBarry Smith 
107947c6ae99SBarry Smith   /* determine the matrix preallocation information */
108047c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
108147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
10821321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10831321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
108447c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
10851321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10861321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
108747c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
10881321219cSEthan Coon 	kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
10891321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
109047c6ae99SBarry Smith 
109147c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
109247c6ae99SBarry Smith 
109347c6ae99SBarry Smith 	cnt  = 0;
109447c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
109547c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
109647c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
109747c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
1098aa219208SBarry Smith 		if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
109947c6ae99SBarry Smith 		  cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
110047c6ae99SBarry Smith 		}
110147c6ae99SBarry Smith 	      }
110247c6ae99SBarry Smith 	    }
110347c6ae99SBarry Smith 	  }
110447c6ae99SBarry Smith 	  rows[l] = l + nc*(slot);
110547c6ae99SBarry Smith 	}
1106784ac674SJed Brown 	ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
110747c6ae99SBarry Smith       }
110847c6ae99SBarry Smith     }
110947c6ae99SBarry Smith   }
1110f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
111147c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
111247c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
111347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1114784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1115784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
111647c6ae99SBarry Smith 
111747c6ae99SBarry Smith   /*
111847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
111947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
112047c6ae99SBarry Smith     PETSc ordering.
112147c6ae99SBarry Smith   */
1122fcfd50ebSBarry Smith   if (!da->prealloc_only) {
112347c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
112447c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
112547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
11261321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
11271321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
112847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
11291321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
11301321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
113147c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
11321321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
11331321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
113447c6ae99SBarry Smith 
113547c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
113647c6ae99SBarry Smith 
113747c6ae99SBarry Smith 	  cnt  = 0;
113847c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
113947c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
114047c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
114147c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
1142aa219208SBarry Smith 		  if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
114347c6ae99SBarry Smith 		    cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
114447c6ae99SBarry Smith 		  }
114547c6ae99SBarry Smith 		}
114647c6ae99SBarry Smith 	      }
114747c6ae99SBarry Smith 	    }
114847c6ae99SBarry Smith 	    rows[l]      = l + nc*(slot);
114947c6ae99SBarry Smith 	  }
115047c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
115147c6ae99SBarry Smith 	}
115247c6ae99SBarry Smith       }
115347c6ae99SBarry Smith     }
115447c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
115547c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
115647c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
115747c6ae99SBarry Smith   }
115847c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
115947c6ae99SBarry Smith   PetscFunctionReturn(0);
116047c6ae99SBarry Smith }
116147c6ae99SBarry Smith 
116247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
116347c6ae99SBarry Smith 
116447c6ae99SBarry Smith #undef __FUNCT__
1165ce308e1dSBarry Smith #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ_Fill"
1166ce308e1dSBarry Smith PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM da,Mat J)
1167ce308e1dSBarry Smith {
1168ce308e1dSBarry Smith   PetscErrorCode         ierr;
1169ce308e1dSBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1170ce308e1dSBarry Smith   PetscInt               xs,nx,i,j,gxs,gnx,row,k,l;
1171ce308e1dSBarry Smith   PetscInt               m,dim,s,*cols = PETSC_NULL,nc,col,cnt,*ocols;
1172ce308e1dSBarry Smith   PetscInt               *ofill = dd->ofill;
1173ce308e1dSBarry Smith   PetscScalar            *values;
1174ce308e1dSBarry Smith   DMDABoundaryType       bx;
1175ce308e1dSBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1176ce308e1dSBarry Smith   PetscMPIInt            rank,size;
1177ce308e1dSBarry Smith 
1178ce308e1dSBarry Smith   PetscFunctionBegin;
1179ce308e1dSBarry Smith   if (dd->bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"With fill provided not implemented with periodic boundary conditions");
1180ce308e1dSBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
1181ce308e1dSBarry Smith   ierr = MPI_Comm_size(((PetscObject)da)->comm,&size);CHKERRQ(ierr);
1182ce308e1dSBarry Smith 
1183ce308e1dSBarry Smith   /*
1184ce308e1dSBarry Smith          nc - number of components per grid point
1185ce308e1dSBarry Smith          col - number of colors needed in one direction for single component problem
1186ce308e1dSBarry Smith 
1187ce308e1dSBarry Smith   */
1188ce308e1dSBarry Smith   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
1189ce308e1dSBarry Smith   col    = 2*s + 1;
1190ce308e1dSBarry Smith 
1191ce308e1dSBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1192ce308e1dSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
1193ce308e1dSBarry Smith 
1194ce308e1dSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
1195ce308e1dSBarry Smith   ierr = PetscMalloc2(nx*nc,PetscInt,&cols,nx*nc,PetscInt,&ocols);CHKERRQ(ierr);
1196ce308e1dSBarry Smith   ierr = PetscMemzero(cols,nx*nc*sizeof(PetscInt));CHKERRQ(ierr);
1197ce308e1dSBarry Smith   ierr = PetscMemzero(ocols,nx*nc*sizeof(PetscInt));CHKERRQ(ierr);
1198ce308e1dSBarry Smith 
1199ce308e1dSBarry Smith   /*
1200ce308e1dSBarry Smith         note should be smaller for first and last process with no periodic
1201ce308e1dSBarry Smith         does not handle dfill
1202ce308e1dSBarry Smith   */
1203ce308e1dSBarry Smith   cnt  = 0;
1204ce308e1dSBarry Smith   /* coupling with process to the left */
1205ce308e1dSBarry Smith   for (i=0; i<s; i++) {
1206ce308e1dSBarry Smith     for (j=0; j<nc; j++){
1207ce308e1dSBarry Smith       ocols[cnt] = ((!rank) ? 0 : (s - i)*(ofill[j+1] - ofill[j]));
1208ce308e1dSBarry Smith       cols[cnt]  = nc + (s + i)*(ofill[j+1] - ofill[j]);
1209ce308e1dSBarry Smith       cnt++;
1210ce308e1dSBarry Smith     }
1211ce308e1dSBarry Smith   }
1212ce308e1dSBarry Smith   for (i=s; i<nx-s; i++) {
1213ce308e1dSBarry Smith     for (j=0; j<nc; j++){
1214ce308e1dSBarry Smith       cols[cnt]  = nc + 2*s*(ofill[j+1] - ofill[j]);
1215ce308e1dSBarry Smith       cnt++;
1216ce308e1dSBarry Smith     }
1217ce308e1dSBarry Smith   }
1218ce308e1dSBarry Smith  /* coupling with process to the right */
1219ce308e1dSBarry Smith   for (i=nx-s; i<nx; i++){
1220ce308e1dSBarry Smith     for (j=0; j<nc; j++){
1221ce308e1dSBarry Smith       ocols[cnt] = ((rank == (size-1)) ? 0 : (i - nx + s + 1)*(ofill[j+1] - ofill[j]));
1222ce308e1dSBarry Smith       cols[cnt]  = nc + (s + nx - i - 1)*(ofill[j+1] - ofill[j]);
1223ce308e1dSBarry Smith       cnt++;
1224ce308e1dSBarry Smith     }
1225ce308e1dSBarry Smith   }
1226ce308e1dSBarry Smith 
1227ce308e1dSBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,cols);CHKERRQ(ierr);
1228ce308e1dSBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,cols,0,ocols);CHKERRQ(ierr);
1229ce308e1dSBarry Smith   ierr = PetscFree2(cols,ocols);CHKERRQ(ierr);
1230ce308e1dSBarry Smith 
1231ce308e1dSBarry Smith   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1232ce308e1dSBarry Smith   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1233ce308e1dSBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1234ce308e1dSBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
1235ce308e1dSBarry Smith 
1236ce308e1dSBarry Smith   /*
1237ce308e1dSBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1238ce308e1dSBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1239ce308e1dSBarry Smith     PETSc ordering.
1240ce308e1dSBarry Smith   */
1241ce308e1dSBarry Smith   if (!da->prealloc_only) {
1242ce308e1dSBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1243ce308e1dSBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1244ce308e1dSBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1245ce308e1dSBarry Smith 
1246ce308e1dSBarry Smith     row  = xs*nc;
1247ce308e1dSBarry Smith     /* coupling with process to the left */
1248ce308e1dSBarry Smith     for (i=xs; i<xs+s; i++) {
1249ce308e1dSBarry Smith       for (j=0; j<nc; j++){
1250ce308e1dSBarry Smith         cnt = 0;
1251ce308e1dSBarry Smith         if (rank) {
1252ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1253ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1254ce308e1dSBarry Smith           }
1255ce308e1dSBarry Smith         }
1256ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1257ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1258ce308e1dSBarry Smith         }
1259ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1260ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1261ce308e1dSBarry Smith         }
1262ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1263ce308e1dSBarry Smith         row++;
1264ce308e1dSBarry Smith       }
1265ce308e1dSBarry Smith     }
1266ce308e1dSBarry Smith     for (i=xs+s; i<xs+nx-s; i++) {
1267ce308e1dSBarry Smith       for (j=0; j<nc; j++){
1268ce308e1dSBarry Smith         cnt = 0;
1269ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1270ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1271ce308e1dSBarry Smith         }
1272ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1273ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1274ce308e1dSBarry Smith         }
1275ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1276ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1277ce308e1dSBarry Smith         }
1278ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1279ce308e1dSBarry Smith         row++;
1280ce308e1dSBarry Smith       }
1281ce308e1dSBarry Smith     }
1282ce308e1dSBarry Smith     /* coupling with process to the right */
1283ce308e1dSBarry Smith     for (i=xs+nx-s; i<xs+nx; i++){
1284ce308e1dSBarry Smith       for (j=0; j<nc; j++){
1285ce308e1dSBarry Smith         cnt = 0;
1286ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1287ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1288ce308e1dSBarry Smith         }
1289ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1290ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1291ce308e1dSBarry Smith         }
1292ce308e1dSBarry Smith         if (rank < size-1) {
1293ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1294ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1295ce308e1dSBarry Smith           }
1296ce308e1dSBarry Smith         }
1297ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1298ce308e1dSBarry Smith         row++;
1299ce308e1dSBarry Smith       }
1300ce308e1dSBarry Smith     }
1301ce308e1dSBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1302ce308e1dSBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1303ce308e1dSBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1304ce308e1dSBarry Smith     ierr = PetscFree(cols);CHKERRQ(ierr);
1305ce308e1dSBarry Smith   }
1306ce308e1dSBarry Smith   PetscFunctionReturn(0);
1307ce308e1dSBarry Smith }
1308ce308e1dSBarry Smith 
1309ce308e1dSBarry Smith /* ---------------------------------------------------------------------------------*/
1310ce308e1dSBarry Smith 
1311ce308e1dSBarry Smith #undef __FUNCT__
1312950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ"
1313950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM da,Mat J)
131447c6ae99SBarry Smith {
131547c6ae99SBarry Smith   PetscErrorCode         ierr;
131647c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
131747c6ae99SBarry Smith   PetscInt               m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l;
131847c6ae99SBarry Smith   PetscInt               istart,iend;
131947c6ae99SBarry Smith   PetscScalar            *values;
13201321219cSEthan Coon   DMDABoundaryType       bx;
132147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
132247c6ae99SBarry Smith 
132347c6ae99SBarry Smith   PetscFunctionBegin;
132447c6ae99SBarry Smith   /*
132547c6ae99SBarry Smith          nc - number of components per grid point
132647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
132747c6ae99SBarry Smith 
132847c6ae99SBarry Smith   */
13291321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
133047c6ae99SBarry Smith   col    = 2*s + 1;
133147c6ae99SBarry Smith 
1332aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1333aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
133447c6ae99SBarry Smith 
1335f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
133647c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
133747c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
133847c6ae99SBarry Smith 
13391411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
13401411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1341784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1342784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
134347c6ae99SBarry Smith 
134447c6ae99SBarry Smith   /*
134547c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
134647c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
134747c6ae99SBarry Smith     PETSc ordering.
134847c6ae99SBarry Smith   */
1349fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1350ce308e1dSBarry Smith     ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
135147c6ae99SBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
135247c6ae99SBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
135347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
135447c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
135547c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
135647c6ae99SBarry Smith       slot   = i - gxs;
135747c6ae99SBarry Smith 
135847c6ae99SBarry Smith       cnt  = 0;
135947c6ae99SBarry Smith       for (l=0; l<nc; l++) {
136047c6ae99SBarry Smith 	for (i1=istart; i1<iend+1; i1++) {
136147c6ae99SBarry Smith 	  cols[cnt++] = l + nc*(slot + i1);
136247c6ae99SBarry Smith 	}
136347c6ae99SBarry Smith 	rows[l]      = l + nc*(slot);
136447c6ae99SBarry Smith       }
136547c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
136647c6ae99SBarry Smith     }
136747c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
136847c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
136947c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
137047c6ae99SBarry Smith     ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
1371ce308e1dSBarry Smith   }
137247c6ae99SBarry Smith   PetscFunctionReturn(0);
137347c6ae99SBarry Smith }
137447c6ae99SBarry Smith 
137547c6ae99SBarry Smith #undef __FUNCT__
1376950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIBAIJ"
1377950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM da,Mat J)
137847c6ae99SBarry Smith {
137947c6ae99SBarry Smith   PetscErrorCode         ierr;
138047c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
138147c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
138247c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
138347c6ae99SBarry Smith   MPI_Comm               comm;
138447c6ae99SBarry Smith   PetscScalar            *values;
13851321219cSEthan Coon   DMDABoundaryType       bx,by;
1386aa219208SBarry Smith   DMDAStencilType        st;
138747c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
138847c6ae99SBarry Smith 
138947c6ae99SBarry Smith   PetscFunctionBegin;
139047c6ae99SBarry Smith   /*
139147c6ae99SBarry Smith      nc - number of components per grid point
139247c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
139347c6ae99SBarry Smith   */
13941321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
139547c6ae99SBarry Smith   col = 2*s + 1;
139647c6ae99SBarry Smith 
1397aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1398aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
139947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
140047c6ae99SBarry Smith 
140147c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
140247c6ae99SBarry Smith 
14031411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
14041411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
140547c6ae99SBarry Smith 
140647c6ae99SBarry Smith   /* determine the matrix preallocation information */
140747c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
140847c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
14091321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14101321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
141147c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
14121321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14131321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
141447c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
141547c6ae99SBarry Smith 
141647c6ae99SBarry Smith       /* Find block columns in block row */
141747c6ae99SBarry Smith       cnt  = 0;
141847c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
141947c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1420aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
142147c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
142247c6ae99SBarry Smith           }
142347c6ae99SBarry Smith         }
142447c6ae99SBarry Smith       }
1425784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
142647c6ae99SBarry Smith     }
142747c6ae99SBarry Smith   }
142847c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
142947c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
143047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
143147c6ae99SBarry Smith 
1432784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1433784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
143447c6ae99SBarry Smith 
143547c6ae99SBarry Smith   /*
143647c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
143747c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
143847c6ae99SBarry Smith     PETSc ordering.
143947c6ae99SBarry Smith   */
1440fcfd50ebSBarry Smith   if (!da->prealloc_only) {
144147c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
144247c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
144347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
14441321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14451321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
144647c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
14471321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14481321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
144947c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
145047c6ae99SBarry Smith 	cnt  = 0;
145147c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
145247c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1453aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
145447c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
145547c6ae99SBarry Smith             }
145647c6ae99SBarry Smith           }
145747c6ae99SBarry Smith         }
145847c6ae99SBarry Smith 	ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
145947c6ae99SBarry Smith       }
146047c6ae99SBarry Smith     }
146147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
146247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
146347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
146447c6ae99SBarry Smith   }
146547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
146647c6ae99SBarry Smith   PetscFunctionReturn(0);
146747c6ae99SBarry Smith }
146847c6ae99SBarry Smith 
146947c6ae99SBarry Smith #undef __FUNCT__
1470950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIBAIJ"
1471950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM da,Mat J)
147247c6ae99SBarry Smith {
147347c6ae99SBarry Smith   PetscErrorCode         ierr;
147447c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
147547c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
147647c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
147747c6ae99SBarry Smith   MPI_Comm               comm;
147847c6ae99SBarry Smith   PetscScalar            *values;
14791321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1480aa219208SBarry Smith   DMDAStencilType        st;
148147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
148247c6ae99SBarry Smith 
148347c6ae99SBarry Smith   PetscFunctionBegin;
148447c6ae99SBarry Smith   /*
148547c6ae99SBarry Smith          nc - number of components per grid point
148647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
148747c6ae99SBarry Smith 
148847c6ae99SBarry Smith   */
14891321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
149047c6ae99SBarry Smith   col    = 2*s + 1;
149147c6ae99SBarry Smith 
1492aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1493aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
149447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
149547c6ae99SBarry Smith 
149647c6ae99SBarry Smith   ierr  = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
149747c6ae99SBarry Smith 
14981411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
14991411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
150047c6ae99SBarry Smith 
150147c6ae99SBarry Smith   /* determine the matrix preallocation information */
150247c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
150347c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
15041321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
15051321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
150647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
15071321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15081321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
150947c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
15101321219cSEthan Coon 	kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15111321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
151247c6ae99SBarry Smith 
151347c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
151447c6ae99SBarry Smith 
151547c6ae99SBarry Smith 	/* Find block columns in block row */
151647c6ae99SBarry Smith 	cnt  = 0;
151747c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
151847c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
151947c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1520aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
152147c6ae99SBarry Smith 		cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
152247c6ae99SBarry Smith 	      }
152347c6ae99SBarry Smith 	    }
152447c6ae99SBarry Smith 	  }
152547c6ae99SBarry Smith 	}
1526784ac674SJed Brown 	ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
152747c6ae99SBarry Smith       }
152847c6ae99SBarry Smith     }
152947c6ae99SBarry Smith   }
153047c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
153147c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
153247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
153347c6ae99SBarry Smith 
1534784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1535784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
153647c6ae99SBarry Smith 
153747c6ae99SBarry Smith   /*
153847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
153947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
154047c6ae99SBarry Smith     PETSc ordering.
154147c6ae99SBarry Smith   */
1542fcfd50ebSBarry Smith   if (!da->prealloc_only) {
154347c6ae99SBarry Smith     ierr  = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
154447c6ae99SBarry Smith     ierr  = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
154547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
15461321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
15471321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
154847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
15491321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15501321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
155147c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
15521321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15531321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
155447c6ae99SBarry Smith 
155547c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
155647c6ae99SBarry Smith 
155747c6ae99SBarry Smith 	  cnt  = 0;
155847c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
155947c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
156047c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1561aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
156247c6ae99SBarry Smith                   cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
156347c6ae99SBarry Smith                 }
156447c6ae99SBarry Smith               }
156547c6ae99SBarry Smith             }
156647c6ae99SBarry Smith           }
156747c6ae99SBarry Smith 	  ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
156847c6ae99SBarry Smith 	}
156947c6ae99SBarry Smith       }
157047c6ae99SBarry Smith     }
157147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
157247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157447c6ae99SBarry Smith   }
157547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
157647c6ae99SBarry Smith   PetscFunctionReturn(0);
157747c6ae99SBarry Smith }
157847c6ae99SBarry Smith 
157947c6ae99SBarry Smith #undef __FUNCT__
158047c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular"
158147c6ae99SBarry Smith /*
158247c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
158347c6ae99SBarry Smith   identify in the local ordering with periodic domain.
158447c6ae99SBarry Smith */
158547c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
158647c6ae99SBarry Smith {
158747c6ae99SBarry Smith   PetscErrorCode ierr;
158847c6ae99SBarry Smith   PetscInt       i,n;
158947c6ae99SBarry Smith 
159047c6ae99SBarry Smith   PetscFunctionBegin;
159147c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr);
159247c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr);
159347c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
159447c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
159547c6ae99SBarry Smith   }
159647c6ae99SBarry Smith   *cnt = n;
159747c6ae99SBarry Smith   PetscFunctionReturn(0);
159847c6ae99SBarry Smith }
159947c6ae99SBarry Smith 
160047c6ae99SBarry Smith #undef __FUNCT__
1601950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPISBAIJ"
1602950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM da,Mat J)
160347c6ae99SBarry Smith {
160447c6ae99SBarry Smith   PetscErrorCode         ierr;
160547c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
160647c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
160747c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
160847c6ae99SBarry Smith   MPI_Comm               comm;
160947c6ae99SBarry Smith   PetscScalar            *values;
16101321219cSEthan Coon   DMDABoundaryType       bx,by;
1611aa219208SBarry Smith   DMDAStencilType        st;
161247c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
161347c6ae99SBarry Smith 
161447c6ae99SBarry Smith   PetscFunctionBegin;
161547c6ae99SBarry Smith   /*
161647c6ae99SBarry Smith      nc - number of components per grid point
161747c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
161847c6ae99SBarry Smith   */
16191321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
162047c6ae99SBarry Smith   col = 2*s + 1;
162147c6ae99SBarry Smith 
1622aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1623aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
162447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
162547c6ae99SBarry Smith 
162647c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
162747c6ae99SBarry Smith 
16281411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
16291411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
163047c6ae99SBarry Smith 
163147c6ae99SBarry Smith   /* determine the matrix preallocation information */
1632eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
163347c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
16341321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16351321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
163647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
16371321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16381321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
163947c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
164047c6ae99SBarry Smith 
164147c6ae99SBarry Smith       /* Find block columns in block row */
164247c6ae99SBarry Smith       cnt  = 0;
164347c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
164447c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1645aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) {
164647c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
164747c6ae99SBarry Smith           }
164847c6ae99SBarry Smith         }
164947c6ae99SBarry Smith       }
165047c6ae99SBarry Smith       ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
165147c6ae99SBarry Smith       ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
165247c6ae99SBarry Smith     }
165347c6ae99SBarry Smith   }
165447c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
165547c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
165647c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
165747c6ae99SBarry Smith 
1658784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1659784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
166047c6ae99SBarry Smith 
166147c6ae99SBarry Smith   /*
166247c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
166347c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
166447c6ae99SBarry Smith     PETSc ordering.
166547c6ae99SBarry Smith   */
1666fcfd50ebSBarry Smith   if (!da->prealloc_only) {
166747c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
166847c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
166947c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
16701321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16711321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
167247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
16731321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16741321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
167547c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
167647c6ae99SBarry Smith 
167747c6ae99SBarry Smith         /* Find block columns in block row */
167847c6ae99SBarry Smith         cnt  = 0;
167947c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
168047c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1681aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) {
168247c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
168347c6ae99SBarry Smith             }
168447c6ae99SBarry Smith           }
168547c6ae99SBarry Smith         }
168647c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
168747c6ae99SBarry Smith 	ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
168847c6ae99SBarry Smith       }
168947c6ae99SBarry Smith     }
169047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
169147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169347c6ae99SBarry Smith   }
169447c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
169547c6ae99SBarry Smith   PetscFunctionReturn(0);
169647c6ae99SBarry Smith }
169747c6ae99SBarry Smith 
169847c6ae99SBarry Smith #undef __FUNCT__
1699950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPISBAIJ"
1700950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM da,Mat J)
170147c6ae99SBarry Smith {
170247c6ae99SBarry Smith   PetscErrorCode         ierr;
170347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
170447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
170547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
170647c6ae99SBarry Smith   MPI_Comm               comm;
170747c6ae99SBarry Smith   PetscScalar            *values;
17081321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1709aa219208SBarry Smith   DMDAStencilType        st;
171047c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
171147c6ae99SBarry Smith 
171247c6ae99SBarry Smith   PetscFunctionBegin;
171347c6ae99SBarry Smith   /*
171447c6ae99SBarry Smith      nc - number of components per grid point
171547c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
171647c6ae99SBarry Smith   */
17171321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
171847c6ae99SBarry Smith   col = 2*s + 1;
171947c6ae99SBarry Smith 
1720aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1721aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
172247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
172347c6ae99SBarry Smith 
172447c6ae99SBarry Smith   /* create the matrix */
172547c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
172647c6ae99SBarry Smith 
17271411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
17281411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
172947c6ae99SBarry Smith 
173047c6ae99SBarry Smith   /* determine the matrix preallocation information */
1731eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
173247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
17331321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
17341321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
173547c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
17361321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
17371321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
173847c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
17391321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
17401321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
174147c6ae99SBarry Smith 
174247c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
174347c6ae99SBarry Smith 
174447c6ae99SBarry Smith 	/* Find block columns in block row */
174547c6ae99SBarry Smith 	cnt  = 0;
174647c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
174747c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
174847c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1749aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
175047c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
175147c6ae99SBarry Smith               }
175247c6ae99SBarry Smith             }
175347c6ae99SBarry Smith           }
175447c6ae99SBarry Smith         }
175547c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
175647c6ae99SBarry Smith         ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
175747c6ae99SBarry Smith       }
175847c6ae99SBarry Smith     }
175947c6ae99SBarry Smith   }
176047c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
176147c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
176247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
176347c6ae99SBarry Smith 
1764784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1765784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
176647c6ae99SBarry Smith 
176747c6ae99SBarry Smith   /*
176847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
176947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
177047c6ae99SBarry Smith     PETSc ordering.
177147c6ae99SBarry Smith   */
1772fcfd50ebSBarry Smith   if (!da->prealloc_only) {
177347c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
177447c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
177547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
17761321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
17771321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
177847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
17791321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
17801321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
178147c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
17821321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
17831321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
178447c6ae99SBarry Smith 
178547c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
178647c6ae99SBarry Smith 
178747c6ae99SBarry Smith 	  cnt  = 0;
178847c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
178947c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
179047c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1791aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
179247c6ae99SBarry Smith 		  cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
179347c6ae99SBarry Smith 		}
179447c6ae99SBarry Smith 	      }
179547c6ae99SBarry Smith 	    }
179647c6ae99SBarry Smith 	  }
179747c6ae99SBarry Smith           ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
179847c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
179947c6ae99SBarry Smith 	}
180047c6ae99SBarry Smith       }
180147c6ae99SBarry Smith     }
180247c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
180347c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
180447c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
180547c6ae99SBarry Smith   }
180647c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
180747c6ae99SBarry Smith   PetscFunctionReturn(0);
180847c6ae99SBarry Smith }
180947c6ae99SBarry Smith 
181047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
181147c6ae99SBarry Smith 
181247c6ae99SBarry Smith #undef __FUNCT__
1813950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ_Fill"
1814950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J)
181547c6ae99SBarry Smith {
181647c6ae99SBarry Smith   PetscErrorCode         ierr;
181747c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
181847c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz;
181947c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
182047c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
182147c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
182247c6ae99SBarry Smith   MPI_Comm               comm;
182347c6ae99SBarry Smith   PetscScalar            *values;
18241321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
182547c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1826aa219208SBarry Smith   DMDAStencilType        st;
182747c6ae99SBarry Smith 
182847c6ae99SBarry Smith   PetscFunctionBegin;
182947c6ae99SBarry Smith   /*
183047c6ae99SBarry Smith          nc - number of components per grid point
183147c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
183247c6ae99SBarry Smith 
183347c6ae99SBarry Smith   */
18341321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
183547c6ae99SBarry Smith   col    = 2*s + 1;
183666a15934SBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
183747c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
183866a15934SBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
183947c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
184066a15934SBarry Smith   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
184147c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
184247c6ae99SBarry Smith 
1843aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1844aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
184547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
184647c6ae99SBarry Smith 
184747c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
18481411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
18491411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
185047c6ae99SBarry Smith 
185147c6ae99SBarry Smith   /* determine the matrix preallocation information */
185247c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
185347c6ae99SBarry Smith 
185447c6ae99SBarry Smith 
185547c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
18561321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
18571321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
185847c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
18591321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
18601321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
186147c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
18621321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
18631321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
186447c6ae99SBarry Smith 
186547c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
186647c6ae99SBarry Smith 
186747c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
186847c6ae99SBarry Smith 	  cnt  = 0;
186947c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
187047c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
187147c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
187247c6ae99SBarry Smith 		if (ii || jj || kk) {
1873aa219208SBarry Smith 		  if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
187447c6ae99SBarry Smith 		    for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
187547c6ae99SBarry Smith 		      cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
187647c6ae99SBarry Smith 		  }
187747c6ae99SBarry Smith 		} else {
187847c6ae99SBarry Smith 		  if (dfill) {
187947c6ae99SBarry Smith 		    for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
188047c6ae99SBarry Smith 		      cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
188147c6ae99SBarry Smith 		  } else {
188247c6ae99SBarry Smith 		    for (ifill_col=0; ifill_col<nc; ifill_col++)
188347c6ae99SBarry Smith 		      cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
188447c6ae99SBarry Smith 		  }
188547c6ae99SBarry Smith 		}
188647c6ae99SBarry Smith 	      }
188747c6ae99SBarry Smith 	    }
188847c6ae99SBarry Smith 	  }
188947c6ae99SBarry Smith 	  row  = l + nc*(slot);
1890784ac674SJed Brown 	  ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
189147c6ae99SBarry Smith 	}
189247c6ae99SBarry Smith       }
189347c6ae99SBarry Smith     }
189447c6ae99SBarry Smith   }
189547c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
189647c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
189747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1898784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1899784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
190047c6ae99SBarry Smith 
190147c6ae99SBarry Smith   /*
190247c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
190347c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
190447c6ae99SBarry Smith     PETSc ordering.
190547c6ae99SBarry Smith   */
1906fcfd50ebSBarry Smith   if (!da->prealloc_only) {
190747c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
190847c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
190947c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
19101321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
19111321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
191247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
19131321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
19141321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
191547c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
19161321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
19171321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
191847c6ae99SBarry Smith 
191947c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
192047c6ae99SBarry Smith 
192147c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
192247c6ae99SBarry Smith 	    cnt  = 0;
192347c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
192447c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
192547c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
192647c6ae99SBarry Smith 		  if (ii || jj || kk) {
1927aa219208SBarry Smith 		    if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
192847c6ae99SBarry Smith 		      for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
192947c6ae99SBarry Smith 			cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
193047c6ae99SBarry Smith 		    }
193147c6ae99SBarry Smith 		  } else {
193247c6ae99SBarry Smith 		    if (dfill) {
193347c6ae99SBarry Smith 		      for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
193447c6ae99SBarry Smith 			cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
193547c6ae99SBarry Smith 		    } else {
193647c6ae99SBarry Smith 		      for (ifill_col=0; ifill_col<nc; ifill_col++)
193747c6ae99SBarry Smith 			cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
193847c6ae99SBarry Smith 		    }
193947c6ae99SBarry Smith 		  }
194047c6ae99SBarry Smith 		}
194147c6ae99SBarry Smith 	      }
194247c6ae99SBarry Smith 	    }
194347c6ae99SBarry Smith 	    row  = l + nc*(slot);
194447c6ae99SBarry Smith 	    ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
194547c6ae99SBarry Smith 	  }
194647c6ae99SBarry Smith 	}
194747c6ae99SBarry Smith       }
194847c6ae99SBarry Smith     }
194947c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
195047c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
195147c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
195247c6ae99SBarry Smith   }
195347c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
195447c6ae99SBarry Smith   PetscFunctionReturn(0);
195547c6ae99SBarry Smith }
1956