xref: /petsc/src/dm/impls/da/fdda.c (revision c688c0463f6c25494d0d243b2a43b2b3454e025c)
147c6ae99SBarry Smith 
2b45d2f2cSJed Brown #include <petsc-private/daimpl.h> /*I      "petscdmda.h"     I*/
3c6db04a5SJed Brown #include <petscmat.h>         /*I      "petscmat.h"    I*/
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"
19aa219208SBarry Smith static PetscErrorCode DMDASetBlockFills_Private(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 */
3647c6ae99SBarry Smith   nz = w + 1;
3747c6ae99SBarry Smith   for (i=0; i<w; i++) {
3847c6ae99SBarry Smith     fill[i] = nz;
3947c6ae99SBarry Smith     for (j=0; j<w; j++) {
4047c6ae99SBarry Smith       if (dfill[w*i+j]) {
4147c6ae99SBarry Smith 	fill[nz] = j;
4247c6ae99SBarry Smith 	nz++;
4347c6ae99SBarry Smith       }
4447c6ae99SBarry Smith     }
4547c6ae99SBarry Smith   }
4647c6ae99SBarry Smith   fill[w] = nz;
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith   *rfill = fill;
4947c6ae99SBarry Smith   PetscFunctionReturn(0);
5047c6ae99SBarry Smith }
5147c6ae99SBarry Smith 
5247c6ae99SBarry Smith #undef __FUNCT__
53aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills"
5447c6ae99SBarry Smith /*@
55aa219208SBarry Smith     DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem
56950540a4SJed Brown     of the matrix returned by DMCreateMatrix().
5747c6ae99SBarry Smith 
58aa219208SBarry Smith     Logically Collective on DMDA
5947c6ae99SBarry Smith 
6047c6ae99SBarry Smith     Input Parameter:
6147c6ae99SBarry Smith +   da - the distributed array
6247c6ae99SBarry Smith .   dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block)
6347c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith 
6647c6ae99SBarry Smith     Level: developer
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith     Notes: This only makes sense when you are doing multicomponent problems but using the
6947c6ae99SBarry Smith        MPIAIJ matrix format
7047c6ae99SBarry Smith 
7147c6ae99SBarry Smith            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
7247c6ae99SBarry Smith        representing coupling and 0 entries for missing coupling. For example
7347c6ae99SBarry Smith $             dfill[9] = {1, 0, 0,
7447c6ae99SBarry Smith $                         1, 1, 0,
7547c6ae99SBarry Smith $                         0, 1, 1}
7647c6ae99SBarry Smith        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
7747c6ae99SBarry Smith        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
7847c6ae99SBarry Smith        diagonal block).
7947c6ae99SBarry Smith 
80aa219208SBarry Smith      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
8147c6ae99SBarry Smith      can be represented in the dfill, ofill format
8247c6ae99SBarry Smith 
8347c6ae99SBarry Smith    Contributed by Glenn Hammond
8447c6ae99SBarry Smith 
858ddb5d8bSBarry Smith .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith @*/
887087cfbeSBarry Smith PetscErrorCode  DMDASetBlockFills(DM da,PetscInt *dfill,PetscInt *ofill)
8947c6ae99SBarry Smith {
9047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
9147c6ae99SBarry Smith   PetscErrorCode ierr;
9247c6ae99SBarry Smith 
9347c6ae99SBarry Smith   PetscFunctionBegin;
94aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr);
95aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr);
9647c6ae99SBarry Smith   PetscFunctionReturn(0);
9747c6ae99SBarry Smith }
9847c6ae99SBarry Smith 
9947c6ae99SBarry Smith 
10047c6ae99SBarry Smith #undef __FUNCT__
101e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA"
10219fd82e9SBarry Smith PetscErrorCode  DMCreateColoring_DA(DM da,ISColoringType ctype,MatType mtype,ISColoring *coloring)
10347c6ae99SBarry Smith {
10447c6ae99SBarry Smith   PetscErrorCode   ierr;
10547c6ae99SBarry Smith   PetscInt         dim,m,n,p,nc;
1061321219cSEthan Coon   DMDABoundaryType bx,by,bz;
10747c6ae99SBarry Smith   MPI_Comm         comm;
10847c6ae99SBarry Smith   PetscMPIInt      size;
10947c6ae99SBarry Smith   PetscBool        isBAIJ;
11047c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
11147c6ae99SBarry Smith 
11247c6ae99SBarry Smith   PetscFunctionBegin;
11347c6ae99SBarry Smith   /*
11447c6ae99SBarry Smith                                   m
11547c6ae99SBarry Smith           ------------------------------------------------------
11647c6ae99SBarry Smith          |                                                     |
11747c6ae99SBarry Smith          |                                                     |
11847c6ae99SBarry Smith          |               ----------------------                |
11947c6ae99SBarry Smith          |               |                    |                |
12047c6ae99SBarry Smith       n  |           yn  |                    |                |
12147c6ae99SBarry Smith          |               |                    |                |
12247c6ae99SBarry Smith          |               .---------------------                |
12347c6ae99SBarry Smith          |             (xs,ys)     xn                          |
12447c6ae99SBarry Smith          |            .                                        |
12547c6ae99SBarry Smith          |         (gxs,gys)                                   |
12647c6ae99SBarry Smith          |                                                     |
12747c6ae99SBarry Smith           -----------------------------------------------------
12847c6ae99SBarry Smith   */
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith   /*
13147c6ae99SBarry Smith          nc - number of components per grid point
13247c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
13347c6ae99SBarry Smith 
13447c6ae99SBarry Smith   */
1351321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&bx,&by,&bz,0);CHKERRQ(ierr);
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
13847c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13947c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED){
14047c6ae99SBarry Smith     if (size == 1) {
14147c6ae99SBarry Smith       ctype = IS_COLORING_GLOBAL;
14247c6ae99SBarry Smith     } else if (dim > 1){
1431321219cSEthan Coon       if ((m==1 && bx == DMDA_BOUNDARY_PERIODIC) || (n==1 && by == DMDA_BOUNDARY_PERIODIC) || (p==1 && bz == DMDA_BOUNDARY_PERIODIC)){
14447c6ae99SBarry 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");
14547c6ae99SBarry Smith       }
14647c6ae99SBarry Smith     }
14747c6ae99SBarry Smith   }
14847c6ae99SBarry Smith 
149aa219208SBarry Smith   /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ
15047c6ae99SBarry Smith      matrices is for the blocks, not the individual matrix elements  */
1514833614aSSatish Balay   ierr = PetscStrcmp(mtype,MATBAIJ,&isBAIJ);CHKERRQ(ierr);
1524833614aSSatish Balay   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);}
15347c6ae99SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);}
15447c6ae99SBarry Smith   if (isBAIJ) {
15547c6ae99SBarry Smith     dd->w = 1;
15647c6ae99SBarry Smith     dd->xs = dd->xs/nc;
15747c6ae99SBarry Smith     dd->xe = dd->xe/nc;
15847c6ae99SBarry Smith     dd->Xs = dd->Xs/nc;
15947c6ae99SBarry Smith     dd->Xe = dd->Xe/nc;
16047c6ae99SBarry Smith   }
16147c6ae99SBarry Smith 
16247c6ae99SBarry Smith   /*
163aa219208SBarry Smith      We do not provide a getcoloring function in the DMDA operations because
164aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
16547c6ae99SBarry Smith    more low-level then matrices.
16647c6ae99SBarry Smith   */
16747c6ae99SBarry Smith   if (dim == 1) {
168e727c939SJed Brown     ierr = DMCreateColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
16947c6ae99SBarry Smith   } else if (dim == 2) {
170e727c939SJed Brown     ierr =  DMCreateColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
17147c6ae99SBarry Smith   } else if (dim == 3) {
172e727c939SJed Brown     ierr =  DMCreateColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
17371cd77b2SBarry 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);
17447c6ae99SBarry Smith   if (isBAIJ) {
17547c6ae99SBarry Smith     dd->w = nc;
17647c6ae99SBarry Smith     dd->xs = dd->xs*nc;
17747c6ae99SBarry Smith     dd->xe = dd->xe*nc;
17847c6ae99SBarry Smith     dd->Xs = dd->Xs*nc;
17947c6ae99SBarry Smith     dd->Xe = dd->Xe*nc;
18047c6ae99SBarry Smith   }
18147c6ae99SBarry Smith   PetscFunctionReturn(0);
18247c6ae99SBarry Smith }
18347c6ae99SBarry Smith 
18447c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
18547c6ae99SBarry Smith 
18647c6ae99SBarry Smith #undef __FUNCT__
187e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_MPIAIJ"
188e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
18947c6ae99SBarry Smith {
19047c6ae99SBarry Smith   PetscErrorCode   ierr;
19147c6ae99SBarry Smith   PetscInt         xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col;
19247c6ae99SBarry Smith   PetscInt         ncolors;
19347c6ae99SBarry Smith   MPI_Comm         comm;
1941321219cSEthan Coon   DMDABoundaryType bx,by;
195aa219208SBarry Smith   DMDAStencilType  st;
19647c6ae99SBarry Smith   ISColoringValue  *colors;
19747c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
19847c6ae99SBarry Smith 
19947c6ae99SBarry Smith   PetscFunctionBegin;
20047c6ae99SBarry Smith   /*
20147c6ae99SBarry Smith          nc - number of components per grid point
20247c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
20347c6ae99SBarry Smith 
20447c6ae99SBarry Smith   */
2051321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
20647c6ae99SBarry Smith   col    = 2*s + 1;
207aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
208aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
20947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
21047c6ae99SBarry Smith 
21147c6ae99SBarry Smith   /* special case as taught to us by Paul Hovland */
212aa219208SBarry Smith   if (st == DMDA_STENCIL_STAR && s == 1) {
213e727c939SJed Brown     ierr = DMCreateColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
21447c6ae99SBarry Smith   } else {
21547c6ae99SBarry Smith 
21666a15934SBarry 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\
21747c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", m, col);
21866a15934SBarry 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\
21947c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", n, col);
22047c6ae99SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
22147c6ae99SBarry Smith       if (!dd->localcoloring) {
22247c6ae99SBarry Smith 	ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
22347c6ae99SBarry Smith 	ii = 0;
22447c6ae99SBarry Smith 	for (j=ys; j<ys+ny; j++) {
22547c6ae99SBarry Smith 	  for (i=xs; i<xs+nx; i++) {
22647c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
22747c6ae99SBarry Smith 	      colors[ii++] = k + nc*((i % col) + col*(j % col));
22847c6ae99SBarry Smith 	    }
22947c6ae99SBarry Smith 	  }
23047c6ae99SBarry Smith 	}
23147c6ae99SBarry Smith         ncolors = nc + nc*(col-1 + col*(col-1));
23247c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
23347c6ae99SBarry Smith       }
23447c6ae99SBarry Smith       *coloring = dd->localcoloring;
23547c6ae99SBarry Smith     } else if (ctype == IS_COLORING_GHOSTED) {
23647c6ae99SBarry Smith       if (!dd->ghostedcoloring) {
23747c6ae99SBarry Smith 	ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
23847c6ae99SBarry Smith 	ii = 0;
23947c6ae99SBarry Smith 	for (j=gys; j<gys+gny; j++) {
24047c6ae99SBarry Smith 	  for (i=gxs; i<gxs+gnx; i++) {
24147c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
24247c6ae99SBarry Smith 	      /* the complicated stuff is to handle periodic boundaries */
24347c6ae99SBarry Smith 	      colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
24447c6ae99SBarry Smith 	    }
24547c6ae99SBarry Smith 	  }
24647c6ae99SBarry Smith 	}
24747c6ae99SBarry Smith         ncolors = nc + nc*(col - 1 + col*(col-1));
24847c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
24947c6ae99SBarry Smith         /* PetscIntView(ncolors,(PetscInt *)colors,0); */
25047c6ae99SBarry Smith 
25147c6ae99SBarry Smith 	ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
25247c6ae99SBarry Smith       }
25347c6ae99SBarry Smith       *coloring = dd->ghostedcoloring;
25447c6ae99SBarry Smith     } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
25547c6ae99SBarry Smith   }
25647c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
25747c6ae99SBarry Smith   PetscFunctionReturn(0);
25847c6ae99SBarry Smith }
25947c6ae99SBarry Smith 
26047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
26147c6ae99SBarry Smith 
26247c6ae99SBarry Smith #undef __FUNCT__
263e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_3d_MPIAIJ"
264e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
26547c6ae99SBarry Smith {
26647c6ae99SBarry Smith   PetscErrorCode    ierr;
26747c6ae99SBarry 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;
26847c6ae99SBarry Smith   PetscInt          ncolors;
26947c6ae99SBarry Smith   MPI_Comm          comm;
2701321219cSEthan Coon   DMDABoundaryType  bx,by,bz;
271aa219208SBarry Smith   DMDAStencilType   st;
27247c6ae99SBarry Smith   ISColoringValue   *colors;
27347c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
27447c6ae99SBarry Smith 
27547c6ae99SBarry Smith   PetscFunctionBegin;
27647c6ae99SBarry Smith   /*
27747c6ae99SBarry Smith          nc - number of components per grid point
27847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
27947c6ae99SBarry Smith 
28047c6ae99SBarry Smith   */
2811321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
28247c6ae99SBarry Smith   col    = 2*s + 1;
28366a15934SBarry 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\
28447c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
28566a15934SBarry 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\
28647c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
28766a15934SBarry 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\
28847c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
28947c6ae99SBarry Smith 
290aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
291aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
29247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
29347c6ae99SBarry Smith 
29447c6ae99SBarry Smith   /* create the coloring */
29547c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
29647c6ae99SBarry Smith     if (!dd->localcoloring) {
29747c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
29847c6ae99SBarry Smith       ii = 0;
29947c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
30047c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
30147c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
30247c6ae99SBarry Smith             for (l=0; l<nc; l++) {
30347c6ae99SBarry Smith               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
30447c6ae99SBarry Smith             }
30547c6ae99SBarry Smith           }
30647c6ae99SBarry Smith         }
30747c6ae99SBarry Smith       }
30847c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
30947c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr);
31047c6ae99SBarry Smith     }
31147c6ae99SBarry Smith     *coloring = dd->localcoloring;
31247c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
31347c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
31447c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
31547c6ae99SBarry Smith       ii = 0;
31647c6ae99SBarry Smith       for (k=gzs; k<gzs+gnz; k++) {
31747c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
31847c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
31947c6ae99SBarry Smith             for (l=0; l<nc; l++) {
32047c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
32147c6ae99SBarry Smith               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
32247c6ae99SBarry Smith             }
32347c6ae99SBarry Smith           }
32447c6ae99SBarry Smith         }
32547c6ae99SBarry Smith       }
32647c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
32747c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
32847c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
32947c6ae99SBarry Smith     }
33047c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
33147c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
33247c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
33347c6ae99SBarry Smith   PetscFunctionReturn(0);
33447c6ae99SBarry Smith }
33547c6ae99SBarry Smith 
33647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
33747c6ae99SBarry Smith 
33847c6ae99SBarry Smith #undef __FUNCT__
339e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_1d_MPIAIJ"
340e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
34147c6ae99SBarry Smith {
34247c6ae99SBarry Smith   PetscErrorCode    ierr;
34347c6ae99SBarry Smith   PetscInt          xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
34447c6ae99SBarry Smith   PetscInt          ncolors;
34547c6ae99SBarry Smith   MPI_Comm          comm;
3461321219cSEthan Coon   DMDABoundaryType  bx;
34747c6ae99SBarry Smith   ISColoringValue   *colors;
34847c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
34947c6ae99SBarry Smith 
35047c6ae99SBarry Smith   PetscFunctionBegin;
35147c6ae99SBarry Smith   /*
35247c6ae99SBarry Smith          nc - number of components per grid point
35347c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
35447c6ae99SBarry Smith 
35547c6ae99SBarry Smith   */
3561321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
35747c6ae99SBarry Smith   col    = 2*s + 1;
35847c6ae99SBarry Smith 
35966a15934SBarry 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\
36031e6f798SBarry Smith                                                           by 2*stencil_width + 1 %d\n",(int)m,(int)col);
36147c6ae99SBarry Smith 
362aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
363aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
36447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
36547c6ae99SBarry Smith 
36647c6ae99SBarry Smith   /* create the coloring */
36747c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
36847c6ae99SBarry Smith     if (!dd->localcoloring) {
36947c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
37047c6ae99SBarry Smith       i1 = 0;
37147c6ae99SBarry Smith       for (i=xs; i<xs+nx; i++) {
37247c6ae99SBarry Smith         for (l=0; l<nc; l++) {
37347c6ae99SBarry Smith           colors[i1++] = l + nc*(i % col);
37447c6ae99SBarry Smith         }
37547c6ae99SBarry Smith       }
37647c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
37747c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr);
37847c6ae99SBarry Smith     }
37947c6ae99SBarry Smith     *coloring = dd->localcoloring;
38047c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
38147c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
38247c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
38347c6ae99SBarry Smith       i1 = 0;
38447c6ae99SBarry Smith       for (i=gxs; i<gxs+gnx; i++) {
38547c6ae99SBarry Smith         for (l=0; l<nc; l++) {
38647c6ae99SBarry Smith           /* the complicated stuff is to handle periodic boundaries */
38747c6ae99SBarry Smith           colors[i1++] = l + nc*(SetInRange(i,m) % col);
38847c6ae99SBarry Smith         }
38947c6ae99SBarry Smith       }
39047c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
39147c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
39247c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
39347c6ae99SBarry Smith     }
39447c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
39547c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
39647c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
39747c6ae99SBarry Smith   PetscFunctionReturn(0);
39847c6ae99SBarry Smith }
39947c6ae99SBarry Smith 
40047c6ae99SBarry Smith #undef __FUNCT__
401e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_5pt_MPIAIJ"
402e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
40347c6ae99SBarry Smith {
40447c6ae99SBarry Smith   PetscErrorCode    ierr;
40547c6ae99SBarry Smith   PetscInt          xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
40647c6ae99SBarry Smith   PetscInt          ncolors;
40747c6ae99SBarry Smith   MPI_Comm          comm;
4081321219cSEthan Coon   DMDABoundaryType  bx,by;
40947c6ae99SBarry Smith   ISColoringValue   *colors;
41047c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
41147c6ae99SBarry Smith 
41247c6ae99SBarry Smith   PetscFunctionBegin;
41347c6ae99SBarry Smith   /*
41447c6ae99SBarry Smith          nc - number of components per grid point
41547c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
41647c6ae99SBarry Smith 
41747c6ae99SBarry Smith   */
4181321219cSEthan Coon   ierr   = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);CHKERRQ(ierr);
419aa219208SBarry Smith   ierr   = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
420aa219208SBarry Smith   ierr   = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
42147c6ae99SBarry Smith   ierr   = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
42247c6ae99SBarry Smith 
42366a15934SBarry 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");
42466a15934SBarry 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");
42547c6ae99SBarry Smith 
42647c6ae99SBarry Smith   /* create the coloring */
42747c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
42847c6ae99SBarry Smith     if (!dd->localcoloring) {
42947c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
43047c6ae99SBarry Smith       ii = 0;
43147c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
43247c6ae99SBarry Smith 	for (i=xs; i<xs+nx; i++) {
43347c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
43447c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*j+i) % 5);
43547c6ae99SBarry Smith 	  }
43647c6ae99SBarry Smith 	}
43747c6ae99SBarry Smith       }
43847c6ae99SBarry Smith       ncolors = 5*nc;
43947c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
44047c6ae99SBarry Smith     }
44147c6ae99SBarry Smith     *coloring = dd->localcoloring;
44247c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
44347c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
44447c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
44547c6ae99SBarry Smith       ii = 0;
44647c6ae99SBarry Smith       for (j=gys; j<gys+gny; j++) {
44747c6ae99SBarry Smith 	for (i=gxs; i<gxs+gnx; i++) {
44847c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
44947c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
45047c6ae99SBarry Smith 	  }
45147c6ae99SBarry Smith 	}
45247c6ae99SBarry Smith       }
45347c6ae99SBarry Smith       ncolors = 5*nc;
45447c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
45547c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
45647c6ae99SBarry Smith     }
45747c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
45847c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
45947c6ae99SBarry Smith   PetscFunctionReturn(0);
46047c6ae99SBarry Smith }
46147c6ae99SBarry Smith 
46247c6ae99SBarry Smith /* =========================================================================== */
463950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM,Mat);
464950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM,Mat);
465950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM,Mat);
466950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM,Mat);
467950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM,Mat);
468950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM,Mat);
469950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM,Mat);
470950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM,Mat);
471950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM,Mat);
47247c6ae99SBarry Smith 
47347c6ae99SBarry Smith #undef __FUNCT__
474*c688c046SMatthew G Knepley #define __FUNCT__ "MatSetupDM"
47547c6ae99SBarry Smith /*@
476*c688c046SMatthew G Knepley    MatSetupDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix
47747c6ae99SBarry Smith 
47847c6ae99SBarry Smith    Logically Collective on Mat
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith    Input Parameters:
48147c6ae99SBarry Smith +  mat - the matrix
48247c6ae99SBarry Smith -  da - the da
48347c6ae99SBarry Smith 
48447c6ae99SBarry Smith    Level: intermediate
48547c6ae99SBarry Smith 
48647c6ae99SBarry Smith @*/
487*c688c046SMatthew G Knepley PetscErrorCode MatSetupDM(Mat mat,DM da)
48847c6ae99SBarry Smith {
48947c6ae99SBarry Smith   PetscErrorCode ierr;
49047c6ae99SBarry Smith 
49147c6ae99SBarry Smith   PetscFunctionBegin;
49247c6ae99SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
49347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
494*c688c046SMatthew G Knepley   ierr = PetscTryMethod(mat,"MatSetupDM_C",(Mat,DM),(mat,da));CHKERRQ(ierr);
49547c6ae99SBarry Smith   PetscFunctionReturn(0);
49647c6ae99SBarry Smith }
49747c6ae99SBarry Smith 
49847c6ae99SBarry Smith EXTERN_C_BEGIN
49947c6ae99SBarry Smith #undef __FUNCT__
50047c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA"
5017087cfbeSBarry Smith PetscErrorCode  MatView_MPI_DA(Mat A,PetscViewer viewer)
50247c6ae99SBarry Smith {
5039a42bb27SBarry Smith   DM             da;
50447c6ae99SBarry Smith   PetscErrorCode ierr;
50547c6ae99SBarry Smith   const char     *prefix;
50647c6ae99SBarry Smith   Mat            Anatural;
50747c6ae99SBarry Smith   AO             ao;
50847c6ae99SBarry Smith   PetscInt       rstart,rend,*petsc,i;
50947c6ae99SBarry Smith   IS             is;
51047c6ae99SBarry Smith   MPI_Comm       comm;
51174388724SJed Brown   PetscViewerFormat format;
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith   PetscFunctionBegin;
51474388724SJed Brown   /* Check whether we are just printing info, in which case MatView() already viewed everything we wanted to view */
51574388724SJed Brown   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
51674388724SJed Brown   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
51774388724SJed Brown 
51847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
519*c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
520aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
52147c6ae99SBarry Smith 
522aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
52347c6ae99SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
52447c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr);
52547c6ae99SBarry Smith   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
52647c6ae99SBarry Smith   ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr);
52747c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
52847c6ae99SBarry Smith 
52947c6ae99SBarry Smith   /* call viewer on natural ordering */
53047c6ae99SBarry Smith   ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr);
531fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
53247c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr);
53347c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr);
53447c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr);
53547c6ae99SBarry Smith   ierr = MatView(Anatural,viewer);CHKERRQ(ierr);
536fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
53747c6ae99SBarry Smith   PetscFunctionReturn(0);
53847c6ae99SBarry Smith }
53947c6ae99SBarry Smith EXTERN_C_END
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith EXTERN_C_BEGIN
54247c6ae99SBarry Smith #undef __FUNCT__
54347c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA"
5447087cfbeSBarry Smith PetscErrorCode  MatLoad_MPI_DA(Mat A,PetscViewer viewer)
54547c6ae99SBarry Smith {
5469a42bb27SBarry Smith   DM             da;
54747c6ae99SBarry Smith   PetscErrorCode ierr;
54847c6ae99SBarry Smith   Mat            Anatural,Aapp;
54947c6ae99SBarry Smith   AO             ao;
55047c6ae99SBarry Smith   PetscInt       rstart,rend,*app,i;
55147c6ae99SBarry Smith   IS             is;
55247c6ae99SBarry Smith   MPI_Comm       comm;
55347c6ae99SBarry Smith 
55447c6ae99SBarry Smith   PetscFunctionBegin;
55547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
556*c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
557aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
55847c6ae99SBarry Smith 
55947c6ae99SBarry Smith   /* Load the matrix in natural ordering */
56047c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,&Anatural);CHKERRQ(ierr);
56147c6ae99SBarry Smith   ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr);
56247c6ae99SBarry Smith   ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
56347c6ae99SBarry Smith   ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr);
56447c6ae99SBarry Smith 
56547c6ae99SBarry Smith   /* Map natural ordering to application ordering and create IS */
566aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
56747c6ae99SBarry Smith   ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr);
56847c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr);
56947c6ae99SBarry Smith   for (i=rstart; i<rend; i++) app[i-rstart] = i;
57047c6ae99SBarry Smith   ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr);
57147c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
57247c6ae99SBarry Smith 
57347c6ae99SBarry Smith   /* Do permutation and replace header */
57447c6ae99SBarry Smith   ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr);
57547c6ae99SBarry Smith   ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr);
576fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
577fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
57847c6ae99SBarry Smith   PetscFunctionReturn(0);
57947c6ae99SBarry Smith }
58047c6ae99SBarry Smith EXTERN_C_END
58147c6ae99SBarry Smith 
58247c6ae99SBarry Smith #undef __FUNCT__
583950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA"
58419fd82e9SBarry Smith PetscErrorCode DMCreateMatrix_DA(DM da, MatType mtype,Mat *J)
58547c6ae99SBarry Smith {
58647c6ae99SBarry Smith   PetscErrorCode ierr;
58747c6ae99SBarry Smith   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
58847c6ae99SBarry Smith   Mat            A;
58947c6ae99SBarry Smith   MPI_Comm       comm;
59019fd82e9SBarry Smith   MatType        Atype;
59137d0c07bSMatthew G Knepley   PetscSection   section, sectionGlobal;
59247c6ae99SBarry Smith   void           (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL;
59347c6ae99SBarry Smith   MatType        ttype[256];
59447c6ae99SBarry Smith   PetscBool      flg;
59547c6ae99SBarry Smith   PetscMPIInt    size;
59647c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
59747c6ae99SBarry Smith 
59847c6ae99SBarry Smith   PetscFunctionBegin;
59947c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
60047c6ae99SBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
60147c6ae99SBarry Smith #endif
6025da5aae0SJed Brown   if (!mtype) mtype = MATAIJ;
60347c6ae99SBarry Smith   ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr);
604aa219208SBarry Smith   ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DMDA options","Mat");CHKERRQ(ierr);
605dd85299cSBarry Smith   ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr);
60647c6ae99SBarry Smith   ierr = PetscOptionsEnd();
60747c6ae99SBarry Smith 
60837d0c07bSMatthew G Knepley   ierr = DMGetDefaultSection(da, &section);CHKERRQ(ierr);
60937d0c07bSMatthew G Knepley   if (section) {
61037d0c07bSMatthew G Knepley     PetscInt  bs = -1;
61137d0c07bSMatthew G Knepley     PetscInt  localSize;
61237d0c07bSMatthew G Knepley     PetscBool isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isSymmetric;
61337d0c07bSMatthew G Knepley 
61437d0c07bSMatthew G Knepley     ierr = DMGetDefaultGlobalSection(da, &sectionGlobal);CHKERRQ(ierr);
61537d0c07bSMatthew G Knepley     ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize);CHKERRQ(ierr);
61637d0c07bSMatthew G Knepley     ierr = MatCreate(((PetscObject) da)->comm, J);CHKERRQ(ierr);
61737d0c07bSMatthew G Knepley     ierr = MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
61837d0c07bSMatthew G Knepley     ierr = MatSetType(*J, mtype);CHKERRQ(ierr);
61937d0c07bSMatthew G Knepley     ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
62037d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSHELL, &isShell);CHKERRQ(ierr);
62137d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATBAIJ, &isBlock);CHKERRQ(ierr);
62237d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSEQBAIJ, &isSeqBlock);CHKERRQ(ierr);
62337d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATMPIBAIJ, &isMPIBlock);CHKERRQ(ierr);
62437d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr);
62537d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr);
62637d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr);
62737d0c07bSMatthew G Knepley     /* Check for symmetric storage */
62837d0c07bSMatthew G Knepley     isSymmetric = (PetscBool) (isSymBlock || isSymSeqBlock || isSymMPIBlock);
62937d0c07bSMatthew G Knepley     if (isSymmetric) {
63037d0c07bSMatthew G Knepley       ierr = MatSetOption(*J, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);
63137d0c07bSMatthew G Knepley     }
63237d0c07bSMatthew G Knepley     if (!isShell) {
6334020307fSSatish Balay       /* PetscBool fillMatrix = (PetscBool) !da->prealloc_only; */
63437d0c07bSMatthew G Knepley       PetscInt *dnz, *onz, *dnzu, *onzu, bsLocal;
63537d0c07bSMatthew G Knepley 
63637d0c07bSMatthew G Knepley       if (bs < 0) {
63737d0c07bSMatthew G Knepley         if (isBlock || isSeqBlock || isMPIBlock || isSymBlock || isSymSeqBlock || isSymMPIBlock) {
63837d0c07bSMatthew G Knepley           PetscInt pStart, pEnd, p, dof;
63937d0c07bSMatthew G Knepley 
64037d0c07bSMatthew G Knepley           ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
64137d0c07bSMatthew G Knepley           for (p = pStart; p < pEnd; ++p) {
64237d0c07bSMatthew G Knepley             ierr = PetscSectionGetDof(sectionGlobal, p, &dof);CHKERRQ(ierr);
64337d0c07bSMatthew G Knepley             if (dof) {
64437d0c07bSMatthew G Knepley               bs = dof;
64537d0c07bSMatthew G Knepley               break;
64637d0c07bSMatthew G Knepley             }
64737d0c07bSMatthew G Knepley           }
64837d0c07bSMatthew G Knepley         } else {
64937d0c07bSMatthew G Knepley           bs = 1;
65037d0c07bSMatthew G Knepley         }
65137d0c07bSMatthew G Knepley         /* Must have same blocksize on all procs (some might have no points) */
65237d0c07bSMatthew G Knepley         bsLocal = bs;
65337d0c07bSMatthew G Knepley         ierr = MPI_Allreduce(&bsLocal, &bs, 1, MPIU_INT, MPI_MAX, ((PetscObject) da)->comm);CHKERRQ(ierr);
65437d0c07bSMatthew G Knepley       }
65537d0c07bSMatthew G Knepley       ierr = PetscMalloc4(localSize/bs, PetscInt, &dnz, localSize/bs, PetscInt, &onz, localSize/bs, PetscInt, &dnzu, localSize/bs, PetscInt, &onzu);CHKERRQ(ierr);
65637d0c07bSMatthew G Knepley       ierr = PetscMemzero(dnz,  localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
65737d0c07bSMatthew G Knepley       ierr = PetscMemzero(onz,  localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
65837d0c07bSMatthew G Knepley       ierr = PetscMemzero(dnzu, localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
65937d0c07bSMatthew G Knepley       ierr = PetscMemzero(onzu, localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
66037d0c07bSMatthew G Knepley       /* ierr = DMComplexPreallocateOperator(dm, bs, section, sectionGlobal, dnz, onz, dnzu, onzu, *J, fillMatrix);CHKERRQ(ierr); */
66137d0c07bSMatthew G Knepley       ierr = PetscFree4(dnz, onz, dnzu, onzu);CHKERRQ(ierr);
66237d0c07bSMatthew G Knepley     }
66337d0c07bSMatthew G Knepley   }
66447c6ae99SBarry Smith   /*
66547c6ae99SBarry Smith                                   m
66647c6ae99SBarry Smith           ------------------------------------------------------
66747c6ae99SBarry Smith          |                                                     |
66847c6ae99SBarry Smith          |                                                     |
66947c6ae99SBarry Smith          |               ----------------------                |
67047c6ae99SBarry Smith          |               |                    |                |
67147c6ae99SBarry Smith       n  |           ny  |                    |                |
67247c6ae99SBarry Smith          |               |                    |                |
67347c6ae99SBarry Smith          |               .---------------------                |
67447c6ae99SBarry Smith          |             (xs,ys)     nx                          |
67547c6ae99SBarry Smith          |            .                                        |
67647c6ae99SBarry Smith          |         (gxs,gys)                                   |
67747c6ae99SBarry Smith          |                                                     |
67847c6ae99SBarry Smith           -----------------------------------------------------
67947c6ae99SBarry Smith   */
68047c6ae99SBarry Smith 
68147c6ae99SBarry Smith   /*
68247c6ae99SBarry Smith          nc - number of components per grid point
68347c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
68447c6ae99SBarry Smith 
68547c6ae99SBarry Smith   */
6861321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
687aa219208SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
68847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
68947c6ae99SBarry Smith   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
69047c6ae99SBarry Smith   ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr);
69119fd82e9SBarry Smith   ierr = MatSetType(A,(MatType)ttype);CHKERRQ(ierr);
69295ee5b0eSBarry Smith   ierr = MatSetDM(A,da);CHKERRQ(ierr);
69347c6ae99SBarry Smith   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
69447c6ae99SBarry Smith   ierr = MatGetType(A,&Atype);CHKERRQ(ierr);
69547c6ae99SBarry Smith   /*
696aa219208SBarry Smith      We do not provide a getmatrix function in the DMDA operations because
697aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
69847c6ae99SBarry Smith    more low-level than matrices. This is kind of cheating but, cause sometimes
699aa219208SBarry Smith    we think of DMDA has higher level than matrices.
70047c6ae99SBarry Smith 
70147c6ae99SBarry Smith      We could switch based on Atype (or mtype), but we do not since the
70247c6ae99SBarry Smith    specialized setting routines depend only the particular preallocation
70347c6ae99SBarry Smith    details of the matrix, not the type itself.
70447c6ae99SBarry Smith   */
70547c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
70647c6ae99SBarry Smith   if (!aij) {
70747c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
70847c6ae99SBarry Smith   }
70947c6ae99SBarry Smith   if (!aij) {
71047c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
71147c6ae99SBarry Smith     if (!baij) {
71247c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
71347c6ae99SBarry Smith     }
71447c6ae99SBarry Smith     if (!baij){
71547c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
71647c6ae99SBarry Smith       if (!sbaij) {
71747c6ae99SBarry Smith         ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
71847c6ae99SBarry Smith       }
71947c6ae99SBarry Smith     }
72047c6ae99SBarry Smith   }
72147c6ae99SBarry Smith   if (aij) {
72247c6ae99SBarry Smith     if (dim == 1) {
723950540a4SJed Brown       ierr = DMCreateMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr);
72447c6ae99SBarry Smith     } else if (dim == 2) {
72547c6ae99SBarry Smith       if (dd->ofill) {
726950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
72747c6ae99SBarry Smith       } else {
728950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr);
72947c6ae99SBarry Smith       }
73047c6ae99SBarry Smith     } else if (dim == 3) {
73147c6ae99SBarry Smith       if (dd->ofill) {
732950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
73347c6ae99SBarry Smith       } else {
734950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr);
73547c6ae99SBarry Smith       }
73647c6ae99SBarry Smith     }
73747c6ae99SBarry Smith   } else if (baij) {
73847c6ae99SBarry Smith     if (dim == 2) {
739950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr);
74047c6ae99SBarry Smith     } else if (dim == 3) {
741950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr);
74266a15934SBarry 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);
74347c6ae99SBarry Smith   } else if (sbaij) {
74447c6ae99SBarry Smith     if (dim == 2) {
745950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr);
74647c6ae99SBarry Smith     } else if (dim == 3) {
747950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr);
74866a15934SBarry 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);
749869776cdSLisandro Dalcin   } else {
750869776cdSLisandro Dalcin     ISLocalToGlobalMapping ltog,ltogb;
751869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
752869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
7532949035bSJed Brown     ierr = MatSetUp(A);CHKERRQ(ierr);
754869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMapping(A,ltog,ltog);CHKERRQ(ierr);
755869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMappingBlock(A,ltogb,ltogb);CHKERRQ(ierr);
75647c6ae99SBarry Smith   }
757aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
75847c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
759*c688c046SMatthew G Knepley   ierr = MatSetDM(A,da);CHKERRQ(ierr);
76047c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
76147c6ae99SBarry Smith   if (size > 1) {
76247c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
76347c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);CHKERRQ(ierr);
76447c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);CHKERRQ(ierr);
76547c6ae99SBarry Smith   }
76647c6ae99SBarry Smith   *J = A;
76747c6ae99SBarry Smith   PetscFunctionReturn(0);
76847c6ae99SBarry Smith }
76947c6ae99SBarry Smith 
77047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
77147c6ae99SBarry Smith #undef __FUNCT__
772950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ"
773950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM da,Mat J)
77447c6ae99SBarry Smith {
77547c6ae99SBarry Smith   PetscErrorCode         ierr;
77647c6ae99SBarry 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;
77747c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
77847c6ae99SBarry Smith   MPI_Comm               comm;
77947c6ae99SBarry Smith   PetscScalar            *values;
7801321219cSEthan Coon   DMDABoundaryType       bx,by;
78147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
782aa219208SBarry Smith   DMDAStencilType        st;
78347c6ae99SBarry Smith 
78447c6ae99SBarry Smith   PetscFunctionBegin;
78547c6ae99SBarry Smith   /*
78647c6ae99SBarry Smith          nc - number of components per grid point
78747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
78847c6ae99SBarry Smith 
78947c6ae99SBarry Smith   */
7901321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
79147c6ae99SBarry Smith   col = 2*s + 1;
792aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
793aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
79447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
79547c6ae99SBarry Smith 
79647c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
7971411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
7981411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
79947c6ae99SBarry Smith 
80047c6ae99SBarry Smith   /* determine the matrix preallocation information */
80147c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
80247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
80347c6ae99SBarry Smith 
8041321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8051321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
80647c6ae99SBarry Smith 
80747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
80847c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
80947c6ae99SBarry Smith 
8101321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8111321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
81247c6ae99SBarry Smith 
81347c6ae99SBarry Smith       cnt  = 0;
81447c6ae99SBarry Smith       for (k=0; k<nc; k++) {
81547c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
81647c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
817aa219208SBarry Smith 	    if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
81847c6ae99SBarry Smith 	      cols[cnt++]  = k + nc*(slot + gnx*l + p);
81947c6ae99SBarry Smith 	    }
82047c6ae99SBarry Smith 	  }
82147c6ae99SBarry Smith 	}
82247c6ae99SBarry Smith 	rows[k] = k + nc*(slot);
82347c6ae99SBarry Smith       }
824784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
82547c6ae99SBarry Smith     }
82647c6ae99SBarry Smith   }
827f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
82847c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
82947c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
83047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
83147c6ae99SBarry Smith 
832784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
833784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
83447c6ae99SBarry Smith 
83547c6ae99SBarry Smith   /*
83647c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
83747c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
83847c6ae99SBarry Smith     PETSc ordering.
83947c6ae99SBarry Smith   */
840fcfd50ebSBarry Smith   if (!da->prealloc_only) {
84147c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
84247c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
84347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
84447c6ae99SBarry Smith 
8451321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8461321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
84947c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
85047c6ae99SBarry Smith 
8511321219cSEthan Coon 	lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8521321219cSEthan Coon 	lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith 	cnt  = 0;
85547c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
85647c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
85747c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
858aa219208SBarry Smith 	      if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
85947c6ae99SBarry Smith 		cols[cnt++]  = k + nc*(slot + gnx*l + p);
86047c6ae99SBarry Smith 	      }
86147c6ae99SBarry Smith 	    }
86247c6ae99SBarry Smith 	  }
86347c6ae99SBarry Smith 	  rows[k]      = k + nc*(slot);
86447c6ae99SBarry Smith 	}
86547c6ae99SBarry Smith 	ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
86647c6ae99SBarry Smith       }
86747c6ae99SBarry Smith     }
86847c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
86947c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
87047c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
87147c6ae99SBarry Smith   }
87247c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
87347c6ae99SBarry Smith   PetscFunctionReturn(0);
87447c6ae99SBarry Smith }
87547c6ae99SBarry Smith 
87647c6ae99SBarry Smith #undef __FUNCT__
877950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ_Fill"
878950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J)
87947c6ae99SBarry Smith {
88047c6ae99SBarry Smith   PetscErrorCode         ierr;
88147c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
88247c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p;
88347c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
88447c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
88547c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
88647c6ae99SBarry Smith   MPI_Comm               comm;
88747c6ae99SBarry Smith   PetscScalar            *values;
8881321219cSEthan Coon   DMDABoundaryType       bx,by;
88947c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
890aa219208SBarry Smith   DMDAStencilType        st;
89147c6ae99SBarry Smith 
89247c6ae99SBarry Smith   PetscFunctionBegin;
89347c6ae99SBarry Smith   /*
89447c6ae99SBarry Smith          nc - number of components per grid point
89547c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
89647c6ae99SBarry Smith 
89747c6ae99SBarry Smith   */
8981321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
89947c6ae99SBarry Smith   col = 2*s + 1;
900aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
901aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
90247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
90347c6ae99SBarry Smith 
90447c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
9051411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
9061411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
90747c6ae99SBarry Smith 
90847c6ae99SBarry Smith   /* determine the matrix preallocation information */
90947c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
91047c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
91147c6ae99SBarry Smith 
9121321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9131321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
91447c6ae99SBarry Smith 
91547c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
91647c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
91747c6ae99SBarry Smith 
9181321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9191321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
92047c6ae99SBarry Smith 
92147c6ae99SBarry Smith       for (k=0; k<nc; k++) {
92247c6ae99SBarry Smith         cnt  = 0;
92347c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
92447c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
92547c6ae99SBarry Smith             if (l || p) {
926aa219208SBarry Smith 	      if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
92747c6ae99SBarry Smith                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
92847c6ae99SBarry Smith 		  cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
92947c6ae99SBarry Smith 	      }
93047c6ae99SBarry Smith             } else {
93147c6ae99SBarry Smith 	      if (dfill) {
93247c6ae99SBarry Smith 		for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
93347c6ae99SBarry Smith 		  cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
93447c6ae99SBarry Smith 	      } else {
93547c6ae99SBarry Smith 		for (ifill_col=0; ifill_col<nc; ifill_col++)
93647c6ae99SBarry Smith 		  cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
93747c6ae99SBarry Smith 	      }
93847c6ae99SBarry Smith             }
93947c6ae99SBarry Smith 	  }
94047c6ae99SBarry Smith 	}
94147c6ae99SBarry Smith 	row = k + nc*(slot);
942784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
94347c6ae99SBarry Smith       }
94447c6ae99SBarry Smith     }
94547c6ae99SBarry Smith   }
94647c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
94747c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
94847c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
949784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
950784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
95147c6ae99SBarry Smith 
95247c6ae99SBarry Smith   /*
95347c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
95447c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
95547c6ae99SBarry Smith     PETSc ordering.
95647c6ae99SBarry Smith   */
957fcfd50ebSBarry Smith   if (!da->prealloc_only) {
95847c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
95947c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
96047c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
96147c6ae99SBarry Smith 
9621321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9631321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
96447c6ae99SBarry Smith 
96547c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
96647c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
96747c6ae99SBarry Smith 
9681321219cSEthan Coon 	lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9691321219cSEthan Coon 	lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
97047c6ae99SBarry Smith 
97147c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
97247c6ae99SBarry Smith 	  cnt  = 0;
97347c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
97447c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
97547c6ae99SBarry Smith 	      if (l || p) {
976aa219208SBarry Smith 		if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
97747c6ae99SBarry Smith 		  for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
97847c6ae99SBarry Smith 		    cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
97947c6ae99SBarry Smith 		}
98047c6ae99SBarry Smith 	      } else {
98147c6ae99SBarry Smith 		if (dfill) {
98247c6ae99SBarry Smith 		  for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
98347c6ae99SBarry Smith 		    cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
98447c6ae99SBarry Smith 		} else {
98547c6ae99SBarry Smith 		  for (ifill_col=0; ifill_col<nc; ifill_col++)
98647c6ae99SBarry Smith 		    cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
98747c6ae99SBarry Smith 		}
98847c6ae99SBarry Smith 	      }
98947c6ae99SBarry Smith 	    }
99047c6ae99SBarry Smith 	  }
99147c6ae99SBarry Smith 	  row  = k + nc*(slot);
99247c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
99347c6ae99SBarry Smith 	}
99447c6ae99SBarry Smith       }
99547c6ae99SBarry Smith     }
99647c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
99747c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
99847c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
99947c6ae99SBarry Smith   }
100047c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
100147c6ae99SBarry Smith   PetscFunctionReturn(0);
100247c6ae99SBarry Smith }
100347c6ae99SBarry Smith 
100447c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
100547c6ae99SBarry Smith 
100647c6ae99SBarry Smith #undef __FUNCT__
1007950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ"
1008950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM da,Mat J)
100947c6ae99SBarry Smith {
101047c6ae99SBarry Smith   PetscErrorCode         ierr;
101147c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
101247c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL;
101347c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
101447c6ae99SBarry Smith   MPI_Comm               comm;
101547c6ae99SBarry Smith   PetscScalar            *values;
10161321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
101747c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1018aa219208SBarry Smith   DMDAStencilType        st;
101947c6ae99SBarry Smith 
102047c6ae99SBarry Smith   PetscFunctionBegin;
102147c6ae99SBarry Smith   /*
102247c6ae99SBarry Smith          nc - number of components per grid point
102347c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
102447c6ae99SBarry Smith 
102547c6ae99SBarry Smith   */
10261321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
102747c6ae99SBarry Smith   col    = 2*s + 1;
102847c6ae99SBarry Smith 
1029aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1030aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
103147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
103247c6ae99SBarry Smith 
103347c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
10341411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
10351411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
103647c6ae99SBarry Smith 
103747c6ae99SBarry Smith   /* determine the matrix preallocation information */
103847c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
103947c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
10401321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10411321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
104247c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
10431321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10441321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
104547c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
10461321219cSEthan Coon 	kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
10471321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
104847c6ae99SBarry Smith 
104947c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
105047c6ae99SBarry Smith 
105147c6ae99SBarry Smith 	cnt  = 0;
105247c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
105347c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
105447c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
105547c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
1056aa219208SBarry Smith 		if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
105747c6ae99SBarry Smith 		  cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
105847c6ae99SBarry Smith 		}
105947c6ae99SBarry Smith 	      }
106047c6ae99SBarry Smith 	    }
106147c6ae99SBarry Smith 	  }
106247c6ae99SBarry Smith 	  rows[l] = l + nc*(slot);
106347c6ae99SBarry Smith 	}
1064784ac674SJed Brown 	ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
106547c6ae99SBarry Smith       }
106647c6ae99SBarry Smith     }
106747c6ae99SBarry Smith   }
1068f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
106947c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
107047c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
107147c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1072784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1073784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
107447c6ae99SBarry Smith 
107547c6ae99SBarry Smith   /*
107647c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
107747c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
107847c6ae99SBarry Smith     PETSc ordering.
107947c6ae99SBarry Smith   */
1080fcfd50ebSBarry Smith   if (!da->prealloc_only) {
108147c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
108247c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
108347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
10841321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10851321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
108647c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
10871321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10881321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
108947c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
10901321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
10911321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
109247c6ae99SBarry Smith 
109347c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
109447c6ae99SBarry Smith 
109547c6ae99SBarry Smith 	  cnt  = 0;
109647c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
109747c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
109847c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
109947c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
1100aa219208SBarry Smith 		  if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
110147c6ae99SBarry Smith 		    cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
110247c6ae99SBarry Smith 		  }
110347c6ae99SBarry Smith 		}
110447c6ae99SBarry Smith 	      }
110547c6ae99SBarry Smith 	    }
110647c6ae99SBarry Smith 	    rows[l]      = l + nc*(slot);
110747c6ae99SBarry Smith 	  }
110847c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
110947c6ae99SBarry Smith 	}
111047c6ae99SBarry Smith       }
111147c6ae99SBarry Smith     }
111247c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
111347c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
111447c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
111547c6ae99SBarry Smith   }
111647c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
111747c6ae99SBarry Smith   PetscFunctionReturn(0);
111847c6ae99SBarry Smith }
111947c6ae99SBarry Smith 
112047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
112147c6ae99SBarry Smith 
112247c6ae99SBarry Smith #undef __FUNCT__
1123950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ"
1124950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM da,Mat J)
112547c6ae99SBarry Smith {
112647c6ae99SBarry Smith   PetscErrorCode         ierr;
112747c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
112847c6ae99SBarry Smith   PetscInt               m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l;
112947c6ae99SBarry Smith   PetscInt               istart,iend;
113047c6ae99SBarry Smith   PetscScalar            *values;
11311321219cSEthan Coon   DMDABoundaryType       bx;
113247c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
113347c6ae99SBarry Smith 
113447c6ae99SBarry Smith   PetscFunctionBegin;
113547c6ae99SBarry Smith   /*
113647c6ae99SBarry Smith          nc - number of components per grid point
113747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
113847c6ae99SBarry Smith 
113947c6ae99SBarry Smith   */
11401321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
114147c6ae99SBarry Smith   col    = 2*s + 1;
114247c6ae99SBarry Smith 
1143aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1144aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
114547c6ae99SBarry Smith 
1146f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
114747c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
114847c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
114947c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
115047c6ae99SBarry Smith 
11511411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
11521411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1153784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1154784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
115547c6ae99SBarry Smith 
115647c6ae99SBarry Smith   /*
115747c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
115847c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
115947c6ae99SBarry Smith     PETSc ordering.
116047c6ae99SBarry Smith   */
1161fcfd50ebSBarry Smith   if (!da->prealloc_only) {
116247c6ae99SBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
116347c6ae99SBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
116447c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
116547c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
116647c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
116747c6ae99SBarry Smith       slot   = i - gxs;
116847c6ae99SBarry Smith 
116947c6ae99SBarry Smith       cnt  = 0;
117047c6ae99SBarry Smith       for (l=0; l<nc; l++) {
117147c6ae99SBarry Smith 	for (i1=istart; i1<iend+1; i1++) {
117247c6ae99SBarry Smith 	  cols[cnt++] = l + nc*(slot + i1);
117347c6ae99SBarry Smith 	}
117447c6ae99SBarry Smith 	rows[l]      = l + nc*(slot);
117547c6ae99SBarry Smith       }
117647c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
117747c6ae99SBarry Smith     }
117847c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
117947c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
118047c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
118147c6ae99SBarry Smith   }
118247c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
118347c6ae99SBarry Smith   PetscFunctionReturn(0);
118447c6ae99SBarry Smith }
118547c6ae99SBarry Smith 
118647c6ae99SBarry Smith #undef __FUNCT__
1187950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIBAIJ"
1188950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM da,Mat J)
118947c6ae99SBarry Smith {
119047c6ae99SBarry Smith   PetscErrorCode         ierr;
119147c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
119247c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
119347c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
119447c6ae99SBarry Smith   MPI_Comm               comm;
119547c6ae99SBarry Smith   PetscScalar            *values;
11961321219cSEthan Coon   DMDABoundaryType       bx,by;
1197aa219208SBarry Smith   DMDAStencilType        st;
119847c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
119947c6ae99SBarry Smith 
120047c6ae99SBarry Smith   PetscFunctionBegin;
120147c6ae99SBarry Smith   /*
120247c6ae99SBarry Smith      nc - number of components per grid point
120347c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
120447c6ae99SBarry Smith   */
12051321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
120647c6ae99SBarry Smith   col = 2*s + 1;
120747c6ae99SBarry Smith 
1208aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1209aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
121047c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
121147c6ae99SBarry Smith 
121247c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
121347c6ae99SBarry Smith 
12141411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
12151411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
121647c6ae99SBarry Smith 
121747c6ae99SBarry Smith   /* determine the matrix preallocation information */
121847c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
121947c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
12201321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
12211321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
122247c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
12231321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
12241321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
122547c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
122647c6ae99SBarry Smith 
122747c6ae99SBarry Smith       /* Find block columns in block row */
122847c6ae99SBarry Smith       cnt  = 0;
122947c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
123047c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1231aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
123247c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
123347c6ae99SBarry Smith           }
123447c6ae99SBarry Smith         }
123547c6ae99SBarry Smith       }
1236784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
123747c6ae99SBarry Smith     }
123847c6ae99SBarry Smith   }
123947c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
124047c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
124147c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
124247c6ae99SBarry Smith 
1243784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1244784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
124547c6ae99SBarry Smith 
124647c6ae99SBarry Smith   /*
124747c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
124847c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
124947c6ae99SBarry Smith     PETSc ordering.
125047c6ae99SBarry Smith   */
1251fcfd50ebSBarry Smith   if (!da->prealloc_only) {
125247c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
125347c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
125447c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
12551321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
12561321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
125747c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
12581321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
12591321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
126047c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
126147c6ae99SBarry Smith 	cnt  = 0;
126247c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
126347c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1264aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
126547c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
126647c6ae99SBarry Smith             }
126747c6ae99SBarry Smith           }
126847c6ae99SBarry Smith         }
126947c6ae99SBarry Smith 	ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
127047c6ae99SBarry Smith       }
127147c6ae99SBarry Smith     }
127247c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
127347c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
127447c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
127547c6ae99SBarry Smith   }
127647c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
127747c6ae99SBarry Smith   PetscFunctionReturn(0);
127847c6ae99SBarry Smith }
127947c6ae99SBarry Smith 
128047c6ae99SBarry Smith #undef __FUNCT__
1281950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIBAIJ"
1282950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM da,Mat J)
128347c6ae99SBarry Smith {
128447c6ae99SBarry Smith   PetscErrorCode         ierr;
128547c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
128647c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
128747c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
128847c6ae99SBarry Smith   MPI_Comm               comm;
128947c6ae99SBarry Smith   PetscScalar            *values;
12901321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1291aa219208SBarry Smith   DMDAStencilType        st;
129247c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
129347c6ae99SBarry Smith 
129447c6ae99SBarry Smith   PetscFunctionBegin;
129547c6ae99SBarry Smith   /*
129647c6ae99SBarry Smith          nc - number of components per grid point
129747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
129847c6ae99SBarry Smith 
129947c6ae99SBarry Smith   */
13001321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
130147c6ae99SBarry Smith   col    = 2*s + 1;
130247c6ae99SBarry Smith 
1303aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1304aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
130547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
130647c6ae99SBarry Smith 
130747c6ae99SBarry Smith   ierr  = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
130847c6ae99SBarry Smith 
13091411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
13101411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
131147c6ae99SBarry Smith 
131247c6ae99SBarry Smith   /* determine the matrix preallocation information */
131347c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
131447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
13151321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
13161321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
131747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
13181321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
13191321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
132047c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
13211321219cSEthan Coon 	kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
13221321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
132347c6ae99SBarry Smith 
132447c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
132547c6ae99SBarry Smith 
132647c6ae99SBarry Smith 	/* Find block columns in block row */
132747c6ae99SBarry Smith 	cnt  = 0;
132847c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
132947c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
133047c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1331aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
133247c6ae99SBarry Smith 		cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
133347c6ae99SBarry Smith 	      }
133447c6ae99SBarry Smith 	    }
133547c6ae99SBarry Smith 	  }
133647c6ae99SBarry Smith 	}
1337784ac674SJed Brown 	ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
133847c6ae99SBarry Smith       }
133947c6ae99SBarry Smith     }
134047c6ae99SBarry Smith   }
134147c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
134247c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
134347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
134447c6ae99SBarry Smith 
1345784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1346784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
134747c6ae99SBarry Smith 
134847c6ae99SBarry Smith   /*
134947c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
135047c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
135147c6ae99SBarry Smith     PETSc ordering.
135247c6ae99SBarry Smith   */
1353fcfd50ebSBarry Smith   if (!da->prealloc_only) {
135447c6ae99SBarry Smith     ierr  = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
135547c6ae99SBarry Smith     ierr  = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
135647c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
13571321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
13581321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
135947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
13601321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
13611321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
136247c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
13631321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
13641321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
136547c6ae99SBarry Smith 
136647c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
136747c6ae99SBarry Smith 
136847c6ae99SBarry Smith 	  cnt  = 0;
136947c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
137047c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
137147c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1372aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
137347c6ae99SBarry Smith                   cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
137447c6ae99SBarry Smith                 }
137547c6ae99SBarry Smith               }
137647c6ae99SBarry Smith             }
137747c6ae99SBarry Smith           }
137847c6ae99SBarry Smith 	  ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
137947c6ae99SBarry Smith 	}
138047c6ae99SBarry Smith       }
138147c6ae99SBarry Smith     }
138247c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
138347c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
138447c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
138547c6ae99SBarry Smith   }
138647c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
138747c6ae99SBarry Smith   PetscFunctionReturn(0);
138847c6ae99SBarry Smith }
138947c6ae99SBarry Smith 
139047c6ae99SBarry Smith #undef __FUNCT__
139147c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular"
139247c6ae99SBarry Smith /*
139347c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
139447c6ae99SBarry Smith   identify in the local ordering with periodic domain.
139547c6ae99SBarry Smith */
139647c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
139747c6ae99SBarry Smith {
139847c6ae99SBarry Smith   PetscErrorCode ierr;
139947c6ae99SBarry Smith   PetscInt       i,n;
140047c6ae99SBarry Smith 
140147c6ae99SBarry Smith   PetscFunctionBegin;
140247c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr);
140347c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr);
140447c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
140547c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
140647c6ae99SBarry Smith   }
140747c6ae99SBarry Smith   *cnt = n;
140847c6ae99SBarry Smith   PetscFunctionReturn(0);
140947c6ae99SBarry Smith }
141047c6ae99SBarry Smith 
141147c6ae99SBarry Smith #undef __FUNCT__
1412950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPISBAIJ"
1413950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM da,Mat J)
141447c6ae99SBarry Smith {
141547c6ae99SBarry Smith   PetscErrorCode         ierr;
141647c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
141747c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
141847c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
141947c6ae99SBarry Smith   MPI_Comm               comm;
142047c6ae99SBarry Smith   PetscScalar            *values;
14211321219cSEthan Coon   DMDABoundaryType       bx,by;
1422aa219208SBarry Smith   DMDAStencilType        st;
142347c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
142447c6ae99SBarry Smith 
142547c6ae99SBarry Smith   PetscFunctionBegin;
142647c6ae99SBarry Smith   /*
142747c6ae99SBarry Smith      nc - number of components per grid point
142847c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
142947c6ae99SBarry Smith   */
14301321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
143147c6ae99SBarry Smith   col = 2*s + 1;
143247c6ae99SBarry Smith 
1433aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1434aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
143547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
143647c6ae99SBarry Smith 
143747c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
143847c6ae99SBarry Smith 
14391411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
14401411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
144147c6ae99SBarry Smith 
144247c6ae99SBarry Smith   /* determine the matrix preallocation information */
1443eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
144447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
14451321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14461321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
144747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
14481321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14491321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
145047c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
145147c6ae99SBarry Smith 
145247c6ae99SBarry Smith       /* Find block columns in block row */
145347c6ae99SBarry Smith       cnt  = 0;
145447c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
145547c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1456aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) {
145747c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
145847c6ae99SBarry Smith           }
145947c6ae99SBarry Smith         }
146047c6ae99SBarry Smith       }
146147c6ae99SBarry Smith       ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
146247c6ae99SBarry Smith       ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
146347c6ae99SBarry Smith     }
146447c6ae99SBarry Smith   }
146547c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
146647c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
146747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
146847c6ae99SBarry Smith 
1469784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1470784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
147147c6ae99SBarry Smith 
147247c6ae99SBarry Smith   /*
147347c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
147447c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
147547c6ae99SBarry Smith     PETSc ordering.
147647c6ae99SBarry Smith   */
1477fcfd50ebSBarry Smith   if (!da->prealloc_only) {
147847c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
147947c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
148047c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
14811321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14821321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
148347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
14841321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14851321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
148647c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
148747c6ae99SBarry Smith 
148847c6ae99SBarry Smith         /* Find block columns in block row */
148947c6ae99SBarry Smith         cnt  = 0;
149047c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
149147c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1492aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) {
149347c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
149447c6ae99SBarry Smith             }
149547c6ae99SBarry Smith           }
149647c6ae99SBarry Smith         }
149747c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
149847c6ae99SBarry Smith 	ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
149947c6ae99SBarry Smith       }
150047c6ae99SBarry Smith     }
150147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
150247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
150447c6ae99SBarry Smith   }
150547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
150647c6ae99SBarry Smith   PetscFunctionReturn(0);
150747c6ae99SBarry Smith }
150847c6ae99SBarry Smith 
150947c6ae99SBarry Smith #undef __FUNCT__
1510950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPISBAIJ"
1511950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM da,Mat J)
151247c6ae99SBarry Smith {
151347c6ae99SBarry Smith   PetscErrorCode         ierr;
151447c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
151547c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
151647c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
151747c6ae99SBarry Smith   MPI_Comm               comm;
151847c6ae99SBarry Smith   PetscScalar            *values;
15191321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1520aa219208SBarry Smith   DMDAStencilType        st;
152147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
152247c6ae99SBarry Smith 
152347c6ae99SBarry Smith   PetscFunctionBegin;
152447c6ae99SBarry Smith   /*
152547c6ae99SBarry Smith      nc - number of components per grid point
152647c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
152747c6ae99SBarry Smith   */
15281321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
152947c6ae99SBarry Smith   col = 2*s + 1;
153047c6ae99SBarry Smith 
1531aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1532aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
153347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
153447c6ae99SBarry Smith 
153547c6ae99SBarry Smith   /* create the matrix */
153647c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
153747c6ae99SBarry Smith 
15381411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
15391411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
154047c6ae99SBarry Smith 
154147c6ae99SBarry Smith   /* determine the matrix preallocation information */
1542eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
154347c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
15441321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
15451321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
154647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
15471321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15481321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
154947c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
15501321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15511321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
155247c6ae99SBarry Smith 
155347c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
155447c6ae99SBarry Smith 
155547c6ae99SBarry Smith 	/* Find block columns in block row */
155647c6ae99SBarry Smith 	cnt  = 0;
155747c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
155847c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
155947c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1560aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
156147c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
156247c6ae99SBarry Smith               }
156347c6ae99SBarry Smith             }
156447c6ae99SBarry Smith           }
156547c6ae99SBarry Smith         }
156647c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
156747c6ae99SBarry Smith         ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
156847c6ae99SBarry Smith       }
156947c6ae99SBarry Smith     }
157047c6ae99SBarry Smith   }
157147c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
157247c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
157347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
157447c6ae99SBarry Smith 
1575784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1576784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
157747c6ae99SBarry Smith 
157847c6ae99SBarry Smith   /*
157947c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
158047c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
158147c6ae99SBarry Smith     PETSc ordering.
158247c6ae99SBarry Smith   */
1583fcfd50ebSBarry Smith   if (!da->prealloc_only) {
158447c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
158547c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
158647c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
15871321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
15881321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
158947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
15901321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15911321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
159247c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
15931321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15941321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
159547c6ae99SBarry Smith 
159647c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
159747c6ae99SBarry Smith 
159847c6ae99SBarry Smith 	  cnt  = 0;
159947c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
160047c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
160147c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1602aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
160347c6ae99SBarry Smith 		  cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
160447c6ae99SBarry Smith 		}
160547c6ae99SBarry Smith 	      }
160647c6ae99SBarry Smith 	    }
160747c6ae99SBarry Smith 	  }
160847c6ae99SBarry Smith           ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
160947c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
161047c6ae99SBarry Smith 	}
161147c6ae99SBarry Smith       }
161247c6ae99SBarry Smith     }
161347c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
161447c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
161547c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
161647c6ae99SBarry Smith   }
161747c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
161847c6ae99SBarry Smith   PetscFunctionReturn(0);
161947c6ae99SBarry Smith }
162047c6ae99SBarry Smith 
162147c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
162247c6ae99SBarry Smith 
162347c6ae99SBarry Smith #undef __FUNCT__
1624950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ_Fill"
1625950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J)
162647c6ae99SBarry Smith {
162747c6ae99SBarry Smith   PetscErrorCode         ierr;
162847c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
162947c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz;
163047c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
163147c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
163247c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
163347c6ae99SBarry Smith   MPI_Comm               comm;
163447c6ae99SBarry Smith   PetscScalar            *values;
16351321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
163647c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1637aa219208SBarry Smith   DMDAStencilType        st;
163847c6ae99SBarry Smith 
163947c6ae99SBarry Smith   PetscFunctionBegin;
164047c6ae99SBarry Smith   /*
164147c6ae99SBarry Smith          nc - number of components per grid point
164247c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
164347c6ae99SBarry Smith 
164447c6ae99SBarry Smith   */
16451321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
164647c6ae99SBarry Smith   col    = 2*s + 1;
164766a15934SBarry 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\
164847c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
164966a15934SBarry 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\
165047c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
165166a15934SBarry 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\
165247c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
165347c6ae99SBarry Smith 
1654aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1655aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
165647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
165747c6ae99SBarry Smith 
165847c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
16591411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
16601411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
166147c6ae99SBarry Smith 
166247c6ae99SBarry Smith   /* determine the matrix preallocation information */
166347c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
166447c6ae99SBarry Smith 
166547c6ae99SBarry Smith 
166647c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
16671321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16681321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
166947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
16701321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16711321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
167247c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
16731321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
16741321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
167547c6ae99SBarry Smith 
167647c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
167747c6ae99SBarry Smith 
167847c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
167947c6ae99SBarry Smith 	  cnt  = 0;
168047c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
168147c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
168247c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
168347c6ae99SBarry Smith 		if (ii || jj || kk) {
1684aa219208SBarry Smith 		  if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
168547c6ae99SBarry Smith 		    for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
168647c6ae99SBarry Smith 		      cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
168747c6ae99SBarry Smith 		  }
168847c6ae99SBarry Smith 		} else {
168947c6ae99SBarry Smith 		  if (dfill) {
169047c6ae99SBarry Smith 		    for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
169147c6ae99SBarry Smith 		      cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
169247c6ae99SBarry Smith 		  } else {
169347c6ae99SBarry Smith 		    for (ifill_col=0; ifill_col<nc; ifill_col++)
169447c6ae99SBarry Smith 		      cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
169547c6ae99SBarry Smith 		  }
169647c6ae99SBarry Smith 		}
169747c6ae99SBarry Smith 	      }
169847c6ae99SBarry Smith 	    }
169947c6ae99SBarry Smith 	  }
170047c6ae99SBarry Smith 	  row  = l + nc*(slot);
1701784ac674SJed Brown 	  ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
170247c6ae99SBarry Smith 	}
170347c6ae99SBarry Smith       }
170447c6ae99SBarry Smith     }
170547c6ae99SBarry Smith   }
170647c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
170747c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
170847c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1709784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1710784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
171147c6ae99SBarry Smith 
171247c6ae99SBarry Smith   /*
171347c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
171447c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
171547c6ae99SBarry Smith     PETSc ordering.
171647c6ae99SBarry Smith   */
1717fcfd50ebSBarry Smith   if (!da->prealloc_only) {
171847c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
171947c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
172047c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
17211321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
17221321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
172347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
17241321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
17251321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
172647c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
17271321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
17281321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
172947c6ae99SBarry Smith 
173047c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
173147c6ae99SBarry Smith 
173247c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
173347c6ae99SBarry Smith 	    cnt  = 0;
173447c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
173547c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
173647c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
173747c6ae99SBarry Smith 		  if (ii || jj || kk) {
1738aa219208SBarry Smith 		    if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
173947c6ae99SBarry Smith 		      for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
174047c6ae99SBarry Smith 			cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
174147c6ae99SBarry Smith 		    }
174247c6ae99SBarry Smith 		  } else {
174347c6ae99SBarry Smith 		    if (dfill) {
174447c6ae99SBarry Smith 		      for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
174547c6ae99SBarry Smith 			cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
174647c6ae99SBarry Smith 		    } else {
174747c6ae99SBarry Smith 		      for (ifill_col=0; ifill_col<nc; ifill_col++)
174847c6ae99SBarry Smith 			cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
174947c6ae99SBarry Smith 		    }
175047c6ae99SBarry Smith 		  }
175147c6ae99SBarry Smith 		}
175247c6ae99SBarry Smith 	      }
175347c6ae99SBarry Smith 	    }
175447c6ae99SBarry Smith 	    row  = l + nc*(slot);
175547c6ae99SBarry Smith 	    ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
175647c6ae99SBarry Smith 	  }
175747c6ae99SBarry Smith 	}
175847c6ae99SBarry Smith       }
175947c6ae99SBarry Smith     }
176047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
176147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
176247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
176347c6ae99SBarry Smith   }
176447c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
176547c6ae99SBarry Smith   PetscFunctionReturn(0);
176647c6ae99SBarry Smith }
1767