xref: /petsc/src/dm/impls/da/fdda.c (revision 3c0c59f39ee3bd561a1b4dcbd0e58366bfa60aba)
147c6ae99SBarry Smith 
2c6db04a5SJed Brown #include <private/daimpl.h> /*I      "petscdmda.h"     I*/
3c6db04a5SJed Brown #include <petscmat.h>         /*I      "petscmat.h"    I*/
4c6db04a5SJed Brown #include <private/matimpl.h>
547c6ae99SBarry Smith 
609573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring *);
709573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring *);
809573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring *);
909573ac7SBarry Smith extern PetscErrorCode DMGetColoring_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
5694013140SBarry Smith     of the matrix returned by DMGetMatrix().
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 
85aa219208SBarry Smith .seealso DMGetMatrix(), DMDASetGetMatrix(), DMDASetMatPreallocateOnly()
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__
10194013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA"
1027087cfbeSBarry Smith PetscErrorCode  DMGetColoring_DA(DM da,ISColoringType ctype,const 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  */
15147c6ae99SBarry Smith   ierr = PetscStrcmp(mtype,MATBAIJ,&isBAIJ);CHKERRQ(ierr);
15247c6ae99SBarry Smith   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) {
16894013140SBarry Smith     ierr = DMGetColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
16947c6ae99SBarry Smith   } else if (dim == 2) {
17094013140SBarry Smith     ierr =  DMGetColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
17147c6ae99SBarry Smith   } else if (dim == 3) {
17294013140SBarry Smith     ierr =  DMGetColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
17347c6ae99SBarry Smith   } else {
17447c6ae99SBarry Smith     SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim);
17547c6ae99SBarry Smith   }
17647c6ae99SBarry Smith   if (isBAIJ) {
17747c6ae99SBarry Smith     dd->w = nc;
17847c6ae99SBarry Smith     dd->xs = dd->xs*nc;
17947c6ae99SBarry Smith     dd->xe = dd->xe*nc;
18047c6ae99SBarry Smith     dd->Xs = dd->Xs*nc;
18147c6ae99SBarry Smith     dd->Xe = dd->Xe*nc;
18247c6ae99SBarry Smith   }
18347c6ae99SBarry Smith   PetscFunctionReturn(0);
18447c6ae99SBarry Smith }
18547c6ae99SBarry Smith 
18647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
18747c6ae99SBarry Smith 
18847c6ae99SBarry Smith #undef __FUNCT__
18994013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_2d_MPIAIJ"
19094013140SBarry Smith PetscErrorCode DMGetColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
19147c6ae99SBarry Smith {
19247c6ae99SBarry Smith   PetscErrorCode         ierr;
19347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col;
19447c6ae99SBarry Smith   PetscInt               ncolors;
19547c6ae99SBarry Smith   MPI_Comm               comm;
1961321219cSEthan Coon   DMDABoundaryType       bx,by;
197aa219208SBarry Smith   DMDAStencilType        st;
19847c6ae99SBarry Smith   ISColoringValue        *colors;
19947c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
20047c6ae99SBarry Smith 
20147c6ae99SBarry Smith   PetscFunctionBegin;
20247c6ae99SBarry Smith   /*
20347c6ae99SBarry Smith          nc - number of components per grid point
20447c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
20547c6ae99SBarry Smith 
20647c6ae99SBarry Smith   */
2071321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
20847c6ae99SBarry Smith   col    = 2*s + 1;
209aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
210aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
21147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
21247c6ae99SBarry Smith 
21347c6ae99SBarry Smith   /* special case as taught to us by Paul Hovland */
214aa219208SBarry Smith   if (st == DMDA_STENCIL_STAR && s == 1) {
21594013140SBarry Smith     ierr = DMGetColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
21647c6ae99SBarry Smith   } else {
21747c6ae99SBarry Smith 
2181321219cSEthan Coon     if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){
21947c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\
22047c6ae99SBarry Smith                  by 2*stencil_width + 1 (%d)\n", m, col);
22147c6ae99SBarry Smith     }
2221321219cSEthan Coon     if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){
22347c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\
22447c6ae99SBarry Smith                  by 2*stencil_width + 1 (%d)\n", n, col);
22547c6ae99SBarry Smith     }
22647c6ae99SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
22747c6ae99SBarry Smith       if (!dd->localcoloring) {
22847c6ae99SBarry Smith 	ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
22947c6ae99SBarry Smith 	ii = 0;
23047c6ae99SBarry Smith 	for (j=ys; j<ys+ny; j++) {
23147c6ae99SBarry Smith 	  for (i=xs; i<xs+nx; i++) {
23247c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
23347c6ae99SBarry Smith 	      colors[ii++] = k + nc*((i % col) + col*(j % col));
23447c6ae99SBarry Smith 	    }
23547c6ae99SBarry Smith 	  }
23647c6ae99SBarry Smith 	}
23747c6ae99SBarry Smith         ncolors = nc + nc*(col-1 + col*(col-1));
23847c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
23947c6ae99SBarry Smith       }
24047c6ae99SBarry Smith       *coloring = dd->localcoloring;
24147c6ae99SBarry Smith     } else if (ctype == IS_COLORING_GHOSTED) {
24247c6ae99SBarry Smith       if (!dd->ghostedcoloring) {
24347c6ae99SBarry Smith 	ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
24447c6ae99SBarry Smith 	ii = 0;
24547c6ae99SBarry Smith 	for (j=gys; j<gys+gny; j++) {
24647c6ae99SBarry Smith 	  for (i=gxs; i<gxs+gnx; i++) {
24747c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
24847c6ae99SBarry Smith 	      /* the complicated stuff is to handle periodic boundaries */
24947c6ae99SBarry Smith 	      colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
25047c6ae99SBarry Smith 	    }
25147c6ae99SBarry Smith 	  }
25247c6ae99SBarry Smith 	}
25347c6ae99SBarry Smith         ncolors = nc + nc*(col - 1 + col*(col-1));
25447c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
25547c6ae99SBarry Smith         /* PetscIntView(ncolors,(PetscInt *)colors,0); */
25647c6ae99SBarry Smith 
25747c6ae99SBarry Smith 	ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
25847c6ae99SBarry Smith       }
25947c6ae99SBarry Smith       *coloring = dd->ghostedcoloring;
26047c6ae99SBarry Smith     } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
26147c6ae99SBarry Smith   }
26247c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
26347c6ae99SBarry Smith   PetscFunctionReturn(0);
26447c6ae99SBarry Smith }
26547c6ae99SBarry Smith 
26647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
26747c6ae99SBarry Smith 
26847c6ae99SBarry Smith #undef __FUNCT__
26994013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_3d_MPIAIJ"
27094013140SBarry Smith PetscErrorCode DMGetColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
27147c6ae99SBarry Smith {
27247c6ae99SBarry Smith   PetscErrorCode    ierr;
27347c6ae99SBarry 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;
27447c6ae99SBarry Smith   PetscInt          ncolors;
27547c6ae99SBarry Smith   MPI_Comm          comm;
2761321219cSEthan Coon   DMDABoundaryType  bx,by,bz;
277aa219208SBarry Smith   DMDAStencilType   st;
27847c6ae99SBarry Smith   ISColoringValue   *colors;
27947c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
28047c6ae99SBarry Smith 
28147c6ae99SBarry Smith   PetscFunctionBegin;
28247c6ae99SBarry Smith   /*
28347c6ae99SBarry Smith          nc - number of components per grid point
28447c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
28547c6ae99SBarry Smith 
28647c6ae99SBarry Smith   */
2871321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
28847c6ae99SBarry Smith   col    = 2*s + 1;
2891321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){
29047c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
29147c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
29247c6ae99SBarry Smith   }
2931321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){
29447c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
29547c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
29647c6ae99SBarry Smith   }
2971321219cSEthan Coon   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)){
29847c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
29947c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
30047c6ae99SBarry Smith   }
30147c6ae99SBarry Smith 
302aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
303aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
30447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
30547c6ae99SBarry Smith 
30647c6ae99SBarry Smith   /* create the coloring */
30747c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
30847c6ae99SBarry Smith     if (!dd->localcoloring) {
30947c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
31047c6ae99SBarry Smith       ii = 0;
31147c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
31247c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
31347c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
31447c6ae99SBarry Smith             for (l=0; l<nc; l++) {
31547c6ae99SBarry Smith               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
31647c6ae99SBarry Smith             }
31747c6ae99SBarry Smith           }
31847c6ae99SBarry Smith         }
31947c6ae99SBarry Smith       }
32047c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
32147c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr);
32247c6ae99SBarry Smith     }
32347c6ae99SBarry Smith     *coloring = dd->localcoloring;
32447c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
32547c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
32647c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
32747c6ae99SBarry Smith       ii = 0;
32847c6ae99SBarry Smith       for (k=gzs; k<gzs+gnz; k++) {
32947c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
33047c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
33147c6ae99SBarry Smith             for (l=0; l<nc; l++) {
33247c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
33347c6ae99SBarry Smith               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
33447c6ae99SBarry Smith             }
33547c6ae99SBarry Smith           }
33647c6ae99SBarry Smith         }
33747c6ae99SBarry Smith       }
33847c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
33947c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
34047c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
34147c6ae99SBarry Smith     }
34247c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
34347c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
34447c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
34547c6ae99SBarry Smith   PetscFunctionReturn(0);
34647c6ae99SBarry Smith }
34747c6ae99SBarry Smith 
34847c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
34947c6ae99SBarry Smith 
35047c6ae99SBarry Smith #undef __FUNCT__
35194013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_1d_MPIAIJ"
35294013140SBarry Smith PetscErrorCode DMGetColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
35347c6ae99SBarry Smith {
35447c6ae99SBarry Smith   PetscErrorCode    ierr;
35547c6ae99SBarry Smith   PetscInt          xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
35647c6ae99SBarry Smith   PetscInt          ncolors;
35747c6ae99SBarry Smith   MPI_Comm          comm;
3581321219cSEthan Coon   DMDABoundaryType  bx;
35947c6ae99SBarry Smith   ISColoringValue   *colors;
36047c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
36147c6ae99SBarry Smith 
36247c6ae99SBarry Smith   PetscFunctionBegin;
36347c6ae99SBarry Smith   /*
36447c6ae99SBarry Smith          nc - number of components per grid point
36547c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith   */
3681321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
36947c6ae99SBarry Smith   col    = 2*s + 1;
37047c6ae99SBarry Smith 
3711321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) {
37247c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points is divisible\n\
37347c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
37447c6ae99SBarry Smith   }
37547c6ae99SBarry Smith 
376aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
377aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
37847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
37947c6ae99SBarry Smith 
38047c6ae99SBarry Smith   /* create the coloring */
38147c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
38247c6ae99SBarry Smith     if (!dd->localcoloring) {
38347c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
38447c6ae99SBarry Smith       i1 = 0;
38547c6ae99SBarry Smith       for (i=xs; i<xs+nx; i++) {
38647c6ae99SBarry Smith         for (l=0; l<nc; l++) {
38747c6ae99SBarry Smith           colors[i1++] = l + nc*(i % col);
38847c6ae99SBarry Smith         }
38947c6ae99SBarry Smith       }
39047c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
39147c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr);
39247c6ae99SBarry Smith     }
39347c6ae99SBarry Smith     *coloring = dd->localcoloring;
39447c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
39547c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
39647c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
39747c6ae99SBarry Smith       i1 = 0;
39847c6ae99SBarry Smith       for (i=gxs; i<gxs+gnx; i++) {
39947c6ae99SBarry Smith         for (l=0; l<nc; l++) {
40047c6ae99SBarry Smith           /* the complicated stuff is to handle periodic boundaries */
40147c6ae99SBarry Smith           colors[i1++] = l + nc*(SetInRange(i,m) % col);
40247c6ae99SBarry Smith         }
40347c6ae99SBarry Smith       }
40447c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
40547c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
40647c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
40747c6ae99SBarry Smith     }
40847c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
40947c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
41047c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
41147c6ae99SBarry Smith   PetscFunctionReturn(0);
41247c6ae99SBarry Smith }
41347c6ae99SBarry Smith 
41447c6ae99SBarry Smith #undef __FUNCT__
41594013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_2d_5pt_MPIAIJ"
41694013140SBarry Smith PetscErrorCode DMGetColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
41747c6ae99SBarry Smith {
41847c6ae99SBarry Smith   PetscErrorCode    ierr;
41947c6ae99SBarry Smith   PetscInt          xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
42047c6ae99SBarry Smith   PetscInt          ncolors;
42147c6ae99SBarry Smith   MPI_Comm          comm;
4221321219cSEthan Coon   DMDABoundaryType  bx,by;
42347c6ae99SBarry Smith   ISColoringValue   *colors;
42447c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
42547c6ae99SBarry Smith 
42647c6ae99SBarry Smith   PetscFunctionBegin;
42747c6ae99SBarry Smith   /*
42847c6ae99SBarry Smith          nc - number of components per grid point
42947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
43047c6ae99SBarry Smith 
43147c6ae99SBarry Smith   */
4321321219cSEthan Coon   ierr   = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);CHKERRQ(ierr);
433aa219208SBarry Smith   ierr   = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
434aa219208SBarry Smith   ierr   = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
43547c6ae99SBarry Smith   ierr   = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
43647c6ae99SBarry Smith 
4371321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && (m % 5)){
43847c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
43947c6ae99SBarry Smith                  by 5\n");
44047c6ae99SBarry Smith   }
4411321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC && (n % 5)){
44247c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
44347c6ae99SBarry Smith                  by 5\n");
44447c6ae99SBarry Smith   }
44547c6ae99SBarry Smith 
44647c6ae99SBarry Smith   /* create the coloring */
44747c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
44847c6ae99SBarry Smith     if (!dd->localcoloring) {
44947c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
45047c6ae99SBarry Smith       ii = 0;
45147c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
45247c6ae99SBarry Smith 	for (i=xs; i<xs+nx; i++) {
45347c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
45447c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*j+i) % 5);
45547c6ae99SBarry Smith 	  }
45647c6ae99SBarry Smith 	}
45747c6ae99SBarry Smith       }
45847c6ae99SBarry Smith       ncolors = 5*nc;
45947c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
46047c6ae99SBarry Smith     }
46147c6ae99SBarry Smith     *coloring = dd->localcoloring;
46247c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
46347c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
46447c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
46547c6ae99SBarry Smith       ii = 0;
46647c6ae99SBarry Smith       for (j=gys; j<gys+gny; j++) {
46747c6ae99SBarry Smith 	for (i=gxs; i<gxs+gnx; i++) {
46847c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
46947c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
47047c6ae99SBarry Smith 	  }
47147c6ae99SBarry Smith 	}
47247c6ae99SBarry Smith       }
47347c6ae99SBarry Smith       ncolors = 5*nc;
47447c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
47547c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
47647c6ae99SBarry Smith     }
47747c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
47847c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
47947c6ae99SBarry Smith   PetscFunctionReturn(0);
48047c6ae99SBarry Smith }
48147c6ae99SBarry Smith 
48247c6ae99SBarry Smith /* =========================================================================== */
48309573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_1d_MPIAIJ(DM,Mat);
48409573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ(DM,Mat);
48509573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ_Fill(DM,Mat);
48609573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ(DM,Mat);
48709573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ_Fill(DM,Mat);
48809573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIBAIJ(DM,Mat);
48909573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIBAIJ(DM,Mat);
49009573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPISBAIJ(DM,Mat);
49109573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPISBAIJ(DM,Mat);
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith #undef __FUNCT__
49495ee5b0eSBarry Smith #define __FUNCT__ "MatSetDM"
49547c6ae99SBarry Smith /*@
49695ee5b0eSBarry Smith    MatSetDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix
49747c6ae99SBarry Smith 
49847c6ae99SBarry Smith    Logically Collective on Mat
49947c6ae99SBarry Smith 
50047c6ae99SBarry Smith    Input Parameters:
50147c6ae99SBarry Smith +  mat - the matrix
50247c6ae99SBarry Smith -  da - the da
50347c6ae99SBarry Smith 
50447c6ae99SBarry Smith    Level: intermediate
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith @*/
50795ee5b0eSBarry Smith PetscErrorCode  MatSetDM(Mat mat,DM da)
50847c6ae99SBarry Smith {
50947c6ae99SBarry Smith   PetscErrorCode ierr;
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith   PetscFunctionBegin;
51247c6ae99SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
51347c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
51495ee5b0eSBarry Smith   ierr = PetscTryMethod(mat,"MatSetDM_C",(Mat,DM),(mat,da));CHKERRQ(ierr);
51547c6ae99SBarry Smith   PetscFunctionReturn(0);
51647c6ae99SBarry Smith }
51747c6ae99SBarry Smith 
51847c6ae99SBarry Smith EXTERN_C_BEGIN
51947c6ae99SBarry Smith #undef __FUNCT__
52047c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA"
5217087cfbeSBarry Smith PetscErrorCode  MatView_MPI_DA(Mat A,PetscViewer viewer)
52247c6ae99SBarry Smith {
5239a42bb27SBarry Smith   DM             da;
52447c6ae99SBarry Smith   PetscErrorCode ierr;
52547c6ae99SBarry Smith   const char     *prefix;
52647c6ae99SBarry Smith   Mat            Anatural;
52747c6ae99SBarry Smith   AO             ao;
52847c6ae99SBarry Smith   PetscInt       rstart,rend,*petsc,i;
52947c6ae99SBarry Smith   IS             is;
53047c6ae99SBarry Smith   MPI_Comm       comm;
53147c6ae99SBarry Smith 
53247c6ae99SBarry Smith   PetscFunctionBegin;
53347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
534*3c0c59f3SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"DM",(PetscObject*)&da);CHKERRQ(ierr);
535aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
53647c6ae99SBarry Smith 
537aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
53847c6ae99SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
53947c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr);
54047c6ae99SBarry Smith   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
54147c6ae99SBarry Smith   ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr);
54247c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
54347c6ae99SBarry Smith 
54447c6ae99SBarry Smith   /* call viewer on natural ordering */
54547c6ae99SBarry Smith   ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr);
546fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
54747c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr);
54847c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr);
54947c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr);
55047c6ae99SBarry Smith   ierr = MatView(Anatural,viewer);CHKERRQ(ierr);
551fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
55247c6ae99SBarry Smith   PetscFunctionReturn(0);
55347c6ae99SBarry Smith }
55447c6ae99SBarry Smith EXTERN_C_END
55547c6ae99SBarry Smith 
55647c6ae99SBarry Smith EXTERN_C_BEGIN
55747c6ae99SBarry Smith #undef __FUNCT__
55847c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA"
5597087cfbeSBarry Smith PetscErrorCode  MatLoad_MPI_DA(Mat A,PetscViewer viewer)
56047c6ae99SBarry Smith {
5619a42bb27SBarry Smith   DM             da;
56247c6ae99SBarry Smith   PetscErrorCode ierr;
56347c6ae99SBarry Smith   Mat            Anatural,Aapp;
56447c6ae99SBarry Smith   AO             ao;
56547c6ae99SBarry Smith   PetscInt       rstart,rend,*app,i;
56647c6ae99SBarry Smith   IS             is;
56747c6ae99SBarry Smith   MPI_Comm       comm;
56847c6ae99SBarry Smith 
56947c6ae99SBarry Smith   PetscFunctionBegin;
57047c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
571*3c0c59f3SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"DM",(PetscObject*)&da);CHKERRQ(ierr);
572aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
57347c6ae99SBarry Smith 
57447c6ae99SBarry Smith   /* Load the matrix in natural ordering */
57547c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,&Anatural);CHKERRQ(ierr);
57647c6ae99SBarry Smith   ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr);
57747c6ae99SBarry Smith   ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
57847c6ae99SBarry Smith   ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr);
57947c6ae99SBarry Smith 
58047c6ae99SBarry Smith   /* Map natural ordering to application ordering and create IS */
581aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
58247c6ae99SBarry Smith   ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr);
58347c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr);
58447c6ae99SBarry Smith   for (i=rstart; i<rend; i++) app[i-rstart] = i;
58547c6ae99SBarry Smith   ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr);
58647c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
58747c6ae99SBarry Smith 
58847c6ae99SBarry Smith   /* Do permutation and replace header */
58947c6ae99SBarry Smith   ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr);
59047c6ae99SBarry Smith   ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr);
591fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
592fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
59347c6ae99SBarry Smith   PetscFunctionReturn(0);
59447c6ae99SBarry Smith }
59547c6ae99SBarry Smith EXTERN_C_END
59647c6ae99SBarry Smith 
59747c6ae99SBarry Smith #undef __FUNCT__
59894013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA"
5997087cfbeSBarry Smith PetscErrorCode  DMGetMatrix_DA(DM da, const MatType mtype,Mat *J)
60047c6ae99SBarry Smith {
60147c6ae99SBarry Smith   PetscErrorCode ierr;
60247c6ae99SBarry Smith   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
60347c6ae99SBarry Smith   Mat            A;
60447c6ae99SBarry Smith   MPI_Comm       comm;
60547c6ae99SBarry Smith   const MatType  Atype;
60647c6ae99SBarry Smith   void           (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL;
60747c6ae99SBarry Smith   MatType        ttype[256];
60847c6ae99SBarry Smith   PetscBool      flg;
60947c6ae99SBarry Smith   PetscMPIInt    size;
61047c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
61147c6ae99SBarry Smith 
61247c6ae99SBarry Smith   PetscFunctionBegin;
61347c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
61447c6ae99SBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
61547c6ae99SBarry Smith #endif
6165da5aae0SJed Brown   if (!mtype) mtype = MATAIJ;
61747c6ae99SBarry Smith   ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr);
618aa219208SBarry Smith   ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DMDA options","Mat");CHKERRQ(ierr);
61947c6ae99SBarry Smith   ierr = PetscOptionsList("-da_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr);
62047c6ae99SBarry Smith   ierr = PetscOptionsEnd();
62147c6ae99SBarry Smith 
62247c6ae99SBarry Smith   /*
62347c6ae99SBarry Smith                                   m
62447c6ae99SBarry Smith           ------------------------------------------------------
62547c6ae99SBarry Smith          |                                                     |
62647c6ae99SBarry Smith          |                                                     |
62747c6ae99SBarry Smith          |               ----------------------                |
62847c6ae99SBarry Smith          |               |                    |                |
62947c6ae99SBarry Smith       n  |           ny  |                    |                |
63047c6ae99SBarry Smith          |               |                    |                |
63147c6ae99SBarry Smith          |               .---------------------                |
63247c6ae99SBarry Smith          |             (xs,ys)     nx                          |
63347c6ae99SBarry Smith          |            .                                        |
63447c6ae99SBarry Smith          |         (gxs,gys)                                   |
63547c6ae99SBarry Smith          |                                                     |
63647c6ae99SBarry Smith           -----------------------------------------------------
63747c6ae99SBarry Smith   */
63847c6ae99SBarry Smith 
63947c6ae99SBarry Smith   /*
64047c6ae99SBarry Smith          nc - number of components per grid point
64147c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
64247c6ae99SBarry Smith 
64347c6ae99SBarry Smith   */
6441321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
645aa219208SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
64647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
64747c6ae99SBarry Smith   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
64847c6ae99SBarry Smith   ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr);
64947c6ae99SBarry Smith   ierr = MatSetType(A,(const MatType)ttype);CHKERRQ(ierr);
65095ee5b0eSBarry Smith   ierr = MatSetDM(A,da);CHKERRQ(ierr);
65147c6ae99SBarry Smith   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
65247c6ae99SBarry Smith   ierr = MatGetType(A,&Atype);CHKERRQ(ierr);
65347c6ae99SBarry Smith   /*
654aa219208SBarry Smith      We do not provide a getmatrix function in the DMDA operations because
655aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
65647c6ae99SBarry Smith    more low-level than matrices. This is kind of cheating but, cause sometimes
657aa219208SBarry Smith    we think of DMDA has higher level than matrices.
65847c6ae99SBarry Smith 
65947c6ae99SBarry Smith      We could switch based on Atype (or mtype), but we do not since the
66047c6ae99SBarry Smith    specialized setting routines depend only the particular preallocation
66147c6ae99SBarry Smith    details of the matrix, not the type itself.
66247c6ae99SBarry Smith   */
66347c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
66447c6ae99SBarry Smith   if (!aij) {
66547c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
66647c6ae99SBarry Smith   }
66747c6ae99SBarry Smith   if (!aij) {
66847c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
66947c6ae99SBarry Smith     if (!baij) {
67047c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
67147c6ae99SBarry Smith     }
67247c6ae99SBarry Smith     if (!baij){
67347c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
67447c6ae99SBarry Smith       if (!sbaij) {
67547c6ae99SBarry Smith         ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
67647c6ae99SBarry Smith       }
67747c6ae99SBarry Smith     }
67847c6ae99SBarry Smith   }
67947c6ae99SBarry Smith   if (aij) {
68047c6ae99SBarry Smith     if (dim == 1) {
68194013140SBarry Smith       ierr = DMGetMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr);
68247c6ae99SBarry Smith     } else if (dim == 2) {
68347c6ae99SBarry Smith       if (dd->ofill) {
68494013140SBarry Smith         ierr = DMGetMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
68547c6ae99SBarry Smith       } else {
68694013140SBarry Smith         ierr = DMGetMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr);
68747c6ae99SBarry Smith       }
68847c6ae99SBarry Smith     } else if (dim == 3) {
68947c6ae99SBarry Smith       if (dd->ofill) {
69094013140SBarry Smith         ierr = DMGetMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
69147c6ae99SBarry Smith       } else {
69294013140SBarry Smith         ierr = DMGetMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr);
69347c6ae99SBarry Smith       }
69447c6ae99SBarry Smith     }
69547c6ae99SBarry Smith   } else if (baij) {
69647c6ae99SBarry Smith     if (dim == 2) {
69794013140SBarry Smith       ierr = DMGetMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr);
69847c6ae99SBarry Smith     } else if (dim == 3) {
69994013140SBarry Smith       ierr = DMGetMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr);
70047c6ae99SBarry Smith     } else {
701b17742caSSean Farley       SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \
702b17742caSSean Farley 	       "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim);
70347c6ae99SBarry Smith     }
70447c6ae99SBarry Smith   } else if (sbaij) {
70547c6ae99SBarry Smith     if (dim == 2) {
70694013140SBarry Smith       ierr = DMGetMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr);
70747c6ae99SBarry Smith     } else if (dim == 3) {
70894013140SBarry Smith       ierr = DMGetMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr);
70947c6ae99SBarry Smith     } else {
710b17742caSSean Farley       SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \
711b17742caSSean Farley 	       "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim);
71247c6ae99SBarry Smith     }
71347c6ae99SBarry Smith   }
714aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
71547c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
716*3c0c59f3SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"DM",(PetscObject)da);CHKERRQ(ierr);
71747c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
71847c6ae99SBarry Smith   if (size > 1) {
71947c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
72047c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);CHKERRQ(ierr);
72147c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);CHKERRQ(ierr);
72247c6ae99SBarry Smith   }
72347c6ae99SBarry Smith   *J = A;
72447c6ae99SBarry Smith   PetscFunctionReturn(0);
72547c6ae99SBarry Smith }
72647c6ae99SBarry Smith 
72747c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
72847c6ae99SBarry Smith #undef __FUNCT__
72994013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIAIJ"
73094013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ(DM da,Mat J)
73147c6ae99SBarry Smith {
73247c6ae99SBarry Smith   PetscErrorCode         ierr;
73347c6ae99SBarry 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;
73447c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
73547c6ae99SBarry Smith   MPI_Comm               comm;
73647c6ae99SBarry Smith   PetscScalar            *values;
7371321219cSEthan Coon   DMDABoundaryType       bx,by;
73847c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
739aa219208SBarry Smith   DMDAStencilType        st;
74047c6ae99SBarry Smith 
74147c6ae99SBarry Smith   PetscFunctionBegin;
74247c6ae99SBarry Smith   /*
74347c6ae99SBarry Smith          nc - number of components per grid point
74447c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
74547c6ae99SBarry Smith 
74647c6ae99SBarry Smith   */
7471321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
74847c6ae99SBarry Smith   col = 2*s + 1;
749aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
750aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
75147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
75247c6ae99SBarry Smith 
75347c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
7541411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
7551411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
75647c6ae99SBarry Smith 
75747c6ae99SBarry Smith   /* determine the matrix preallocation information */
75847c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
75947c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
76047c6ae99SBarry Smith 
7611321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
7621321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
76347c6ae99SBarry Smith 
76447c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
76547c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
76647c6ae99SBarry Smith 
7671321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
7681321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
76947c6ae99SBarry Smith 
77047c6ae99SBarry Smith       cnt  = 0;
77147c6ae99SBarry Smith       for (k=0; k<nc; k++) {
77247c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
77347c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
774aa219208SBarry Smith 	    if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
77547c6ae99SBarry Smith 	      cols[cnt++]  = k + nc*(slot + gnx*l + p);
77647c6ae99SBarry Smith 	    }
77747c6ae99SBarry Smith 	  }
77847c6ae99SBarry Smith 	}
77947c6ae99SBarry Smith 	rows[k] = k + nc*(slot);
78047c6ae99SBarry Smith       }
781784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
78247c6ae99SBarry Smith     }
78347c6ae99SBarry Smith   }
78447c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
78547c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
78647c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
78747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
78847c6ae99SBarry Smith 
789784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
790784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith   /*
79347c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
79447c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
79547c6ae99SBarry Smith     PETSc ordering.
79647c6ae99SBarry Smith   */
797fcfd50ebSBarry Smith   if (!da->prealloc_only) {
79847c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
79947c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
80047c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
80147c6ae99SBarry Smith 
8021321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8031321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
80447c6ae99SBarry Smith 
80547c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
80647c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
80747c6ae99SBarry Smith 
8081321219cSEthan Coon 	lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8091321219cSEthan Coon 	lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
81047c6ae99SBarry Smith 
81147c6ae99SBarry Smith 	cnt  = 0;
81247c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
81347c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
81447c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
815aa219208SBarry Smith 	      if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
81647c6ae99SBarry Smith 		cols[cnt++]  = k + nc*(slot + gnx*l + p);
81747c6ae99SBarry Smith 	      }
81847c6ae99SBarry Smith 	    }
81947c6ae99SBarry Smith 	  }
82047c6ae99SBarry Smith 	  rows[k]      = k + nc*(slot);
82147c6ae99SBarry Smith 	}
82247c6ae99SBarry Smith 	ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
82347c6ae99SBarry Smith       }
82447c6ae99SBarry Smith     }
82547c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
82647c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
82747c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
82847c6ae99SBarry Smith   }
82947c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
83047c6ae99SBarry Smith   PetscFunctionReturn(0);
83147c6ae99SBarry Smith }
83247c6ae99SBarry Smith 
83347c6ae99SBarry Smith #undef __FUNCT__
83494013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIAIJ_Fill"
83594013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J)
83647c6ae99SBarry Smith {
83747c6ae99SBarry Smith   PetscErrorCode         ierr;
83847c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
83947c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p;
84047c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
84147c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
84247c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
84347c6ae99SBarry Smith   MPI_Comm               comm;
84447c6ae99SBarry Smith   PetscScalar            *values;
8451321219cSEthan Coon   DMDABoundaryType       bx,by;
84647c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
847aa219208SBarry Smith   DMDAStencilType        st;
84847c6ae99SBarry Smith 
84947c6ae99SBarry Smith   PetscFunctionBegin;
85047c6ae99SBarry Smith   /*
85147c6ae99SBarry Smith          nc - number of components per grid point
85247c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith   */
8551321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
85647c6ae99SBarry Smith   col = 2*s + 1;
857aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
858aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
85947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
8621411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
8631411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
86447c6ae99SBarry Smith 
86547c6ae99SBarry Smith   /* determine the matrix preallocation information */
86647c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
86747c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
86847c6ae99SBarry Smith 
8691321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8701321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
87147c6ae99SBarry Smith 
87247c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
87347c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
87447c6ae99SBarry Smith 
8751321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8761321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
87747c6ae99SBarry Smith 
87847c6ae99SBarry Smith       for (k=0; k<nc; k++) {
87947c6ae99SBarry Smith         cnt  = 0;
88047c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
88147c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
88247c6ae99SBarry Smith             if (l || p) {
883aa219208SBarry Smith 	      if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
88447c6ae99SBarry Smith                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
88547c6ae99SBarry Smith 		  cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
88647c6ae99SBarry Smith 	      }
88747c6ae99SBarry Smith             } else {
88847c6ae99SBarry Smith 	      if (dfill) {
88947c6ae99SBarry Smith 		for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
89047c6ae99SBarry Smith 		  cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
89147c6ae99SBarry Smith 	      } else {
89247c6ae99SBarry Smith 		for (ifill_col=0; ifill_col<nc; ifill_col++)
89347c6ae99SBarry Smith 		  cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
89447c6ae99SBarry Smith 	      }
89547c6ae99SBarry Smith             }
89647c6ae99SBarry Smith 	  }
89747c6ae99SBarry Smith 	}
89847c6ae99SBarry Smith 	row = k + nc*(slot);
899784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
90047c6ae99SBarry Smith       }
90147c6ae99SBarry Smith     }
90247c6ae99SBarry Smith   }
90347c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
90447c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
90547c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
906784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
907784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
90847c6ae99SBarry Smith 
90947c6ae99SBarry Smith   /*
91047c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
91147c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
91247c6ae99SBarry Smith     PETSc ordering.
91347c6ae99SBarry Smith   */
914fcfd50ebSBarry Smith   if (!da->prealloc_only) {
91547c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
91647c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
91747c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
91847c6ae99SBarry Smith 
9191321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9201321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
92147c6ae99SBarry Smith 
92247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
92347c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
92447c6ae99SBarry Smith 
9251321219cSEthan Coon 	lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9261321219cSEthan Coon 	lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
92747c6ae99SBarry Smith 
92847c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
92947c6ae99SBarry Smith 	  cnt  = 0;
93047c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
93147c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
93247c6ae99SBarry Smith 	      if (l || p) {
933aa219208SBarry Smith 		if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
93447c6ae99SBarry Smith 		  for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
93547c6ae99SBarry Smith 		    cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
93647c6ae99SBarry Smith 		}
93747c6ae99SBarry Smith 	      } else {
93847c6ae99SBarry Smith 		if (dfill) {
93947c6ae99SBarry Smith 		  for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
94047c6ae99SBarry Smith 		    cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
94147c6ae99SBarry Smith 		} else {
94247c6ae99SBarry Smith 		  for (ifill_col=0; ifill_col<nc; ifill_col++)
94347c6ae99SBarry Smith 		    cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
94447c6ae99SBarry Smith 		}
94547c6ae99SBarry Smith 	      }
94647c6ae99SBarry Smith 	    }
94747c6ae99SBarry Smith 	  }
94847c6ae99SBarry Smith 	  row  = k + nc*(slot);
94947c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
95047c6ae99SBarry Smith 	}
95147c6ae99SBarry Smith       }
95247c6ae99SBarry Smith     }
95347c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
95447c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
95547c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
95647c6ae99SBarry Smith   }
95747c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
95847c6ae99SBarry Smith   PetscFunctionReturn(0);
95947c6ae99SBarry Smith }
96047c6ae99SBarry Smith 
96147c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
96247c6ae99SBarry Smith 
96347c6ae99SBarry Smith #undef __FUNCT__
96494013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIAIJ"
96594013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ(DM da,Mat J)
96647c6ae99SBarry Smith {
96747c6ae99SBarry Smith   PetscErrorCode         ierr;
96847c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
96947c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL;
97047c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
97147c6ae99SBarry Smith   MPI_Comm               comm;
97247c6ae99SBarry Smith   PetscScalar            *values;
9731321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
97447c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
975aa219208SBarry Smith   DMDAStencilType        st;
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith   PetscFunctionBegin;
97847c6ae99SBarry Smith   /*
97947c6ae99SBarry Smith          nc - number of components per grid point
98047c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
98147c6ae99SBarry Smith 
98247c6ae99SBarry Smith   */
9831321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
98447c6ae99SBarry Smith   col    = 2*s + 1;
98547c6ae99SBarry Smith 
986aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
987aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
98847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
98947c6ae99SBarry Smith 
99047c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
9911411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
9921411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
99347c6ae99SBarry Smith 
99447c6ae99SBarry Smith   /* determine the matrix preallocation information */
99547c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
99647c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
9971321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9981321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
99947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
10001321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10011321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
100247c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
10031321219cSEthan Coon 	kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
10041321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
100547c6ae99SBarry Smith 
100647c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
100747c6ae99SBarry Smith 
100847c6ae99SBarry Smith 	cnt  = 0;
100947c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
101047c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
101147c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
101247c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
1013aa219208SBarry Smith 		if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
101447c6ae99SBarry Smith 		  cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
101547c6ae99SBarry Smith 		}
101647c6ae99SBarry Smith 	      }
101747c6ae99SBarry Smith 	    }
101847c6ae99SBarry Smith 	  }
101947c6ae99SBarry Smith 	  rows[l] = l + nc*(slot);
102047c6ae99SBarry Smith 	}
1021784ac674SJed Brown 	ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
102247c6ae99SBarry Smith       }
102347c6ae99SBarry Smith     }
102447c6ae99SBarry Smith   }
102547c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
102647c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
102747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
102847c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
1029784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1030784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
103147c6ae99SBarry Smith 
103247c6ae99SBarry Smith   /*
103347c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
103447c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
103547c6ae99SBarry Smith     PETSc ordering.
103647c6ae99SBarry Smith   */
1037fcfd50ebSBarry Smith   if (!da->prealloc_only) {
103847c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
103947c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
104047c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
10411321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10421321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
104347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
10441321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10451321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
104647c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
10471321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
10481321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
104947c6ae99SBarry Smith 
105047c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
105147c6ae99SBarry Smith 
105247c6ae99SBarry Smith 	  cnt  = 0;
105347c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
105447c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
105547c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
105647c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
1057aa219208SBarry Smith 		  if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
105847c6ae99SBarry Smith 		    cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
105947c6ae99SBarry Smith 		  }
106047c6ae99SBarry Smith 		}
106147c6ae99SBarry Smith 	      }
106247c6ae99SBarry Smith 	    }
106347c6ae99SBarry Smith 	    rows[l]      = l + nc*(slot);
106447c6ae99SBarry Smith 	  }
106547c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
106647c6ae99SBarry Smith 	}
106747c6ae99SBarry Smith       }
106847c6ae99SBarry Smith     }
106947c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
107047c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
107147c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
107247c6ae99SBarry Smith   }
107347c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
107447c6ae99SBarry Smith   PetscFunctionReturn(0);
107547c6ae99SBarry Smith }
107647c6ae99SBarry Smith 
107747c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
107847c6ae99SBarry Smith 
107947c6ae99SBarry Smith #undef __FUNCT__
108094013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_1d_MPIAIJ"
108194013140SBarry Smith PetscErrorCode DMGetMatrix_DA_1d_MPIAIJ(DM da,Mat J)
108247c6ae99SBarry Smith {
108347c6ae99SBarry Smith   PetscErrorCode         ierr;
108447c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
108547c6ae99SBarry Smith   PetscInt               m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l;
108647c6ae99SBarry Smith   PetscInt               istart,iend;
108747c6ae99SBarry Smith   PetscScalar            *values;
10881321219cSEthan Coon   DMDABoundaryType       bx;
108947c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
109047c6ae99SBarry Smith 
109147c6ae99SBarry Smith   PetscFunctionBegin;
109247c6ae99SBarry Smith   /*
109347c6ae99SBarry Smith          nc - number of components per grid point
109447c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
109547c6ae99SBarry Smith 
109647c6ae99SBarry Smith   */
10971321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
109847c6ae99SBarry Smith   col    = 2*s + 1;
109947c6ae99SBarry Smith 
1100aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1101aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
110247c6ae99SBarry Smith 
110347c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
110447c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
110547c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
110647c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
110747c6ae99SBarry Smith 
11081411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
11091411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1110784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1111784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
111247c6ae99SBarry Smith 
111347c6ae99SBarry Smith   /*
111447c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
111547c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
111647c6ae99SBarry Smith     PETSc ordering.
111747c6ae99SBarry Smith   */
1118fcfd50ebSBarry Smith   if (!da->prealloc_only) {
111947c6ae99SBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
112047c6ae99SBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
112147c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
112247c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
112347c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
112447c6ae99SBarry Smith       slot   = i - gxs;
112547c6ae99SBarry Smith 
112647c6ae99SBarry Smith       cnt  = 0;
112747c6ae99SBarry Smith       for (l=0; l<nc; l++) {
112847c6ae99SBarry Smith 	for (i1=istart; i1<iend+1; i1++) {
112947c6ae99SBarry Smith 	  cols[cnt++] = l + nc*(slot + i1);
113047c6ae99SBarry Smith 	}
113147c6ae99SBarry Smith 	rows[l]      = l + nc*(slot);
113247c6ae99SBarry Smith       }
113347c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
113447c6ae99SBarry Smith     }
113547c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
113647c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
113747c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
113847c6ae99SBarry Smith   }
113947c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
114047c6ae99SBarry Smith   PetscFunctionReturn(0);
114147c6ae99SBarry Smith }
114247c6ae99SBarry Smith 
114347c6ae99SBarry Smith #undef __FUNCT__
114494013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIBAIJ"
114594013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIBAIJ(DM da,Mat J)
114647c6ae99SBarry Smith {
114747c6ae99SBarry Smith   PetscErrorCode         ierr;
114847c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
114947c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
115047c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
115147c6ae99SBarry Smith   MPI_Comm               comm;
115247c6ae99SBarry Smith   PetscScalar            *values;
11531321219cSEthan Coon   DMDABoundaryType       bx,by;
1154aa219208SBarry Smith   DMDAStencilType        st;
115547c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
115647c6ae99SBarry Smith 
115747c6ae99SBarry Smith   PetscFunctionBegin;
115847c6ae99SBarry Smith   /*
115947c6ae99SBarry Smith      nc - number of components per grid point
116047c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
116147c6ae99SBarry Smith   */
11621321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
116347c6ae99SBarry Smith   col = 2*s + 1;
116447c6ae99SBarry Smith 
1165aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1166aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
116747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
116847c6ae99SBarry Smith 
116947c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
117047c6ae99SBarry Smith 
11711411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
11721411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
117347c6ae99SBarry Smith 
117447c6ae99SBarry Smith   /* determine the matrix preallocation information */
117547c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
117647c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
11771321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
11781321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
117947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
11801321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
11811321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
118247c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
118347c6ae99SBarry Smith 
118447c6ae99SBarry Smith       /* Find block columns in block row */
118547c6ae99SBarry Smith       cnt  = 0;
118647c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
118747c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1188aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
118947c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
119047c6ae99SBarry Smith           }
119147c6ae99SBarry Smith         }
119247c6ae99SBarry Smith       }
1193784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
119447c6ae99SBarry Smith     }
119547c6ae99SBarry Smith   }
119647c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
119747c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
119847c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
119947c6ae99SBarry Smith 
1200784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1201784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
120247c6ae99SBarry Smith 
120347c6ae99SBarry Smith   /*
120447c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
120547c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
120647c6ae99SBarry Smith     PETSc ordering.
120747c6ae99SBarry Smith   */
1208fcfd50ebSBarry Smith   if (!da->prealloc_only) {
120947c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
121047c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
121147c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
12121321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
12131321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
121447c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
12151321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
12161321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
121747c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
121847c6ae99SBarry Smith 	cnt  = 0;
121947c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
122047c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1221aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
122247c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
122347c6ae99SBarry Smith             }
122447c6ae99SBarry Smith           }
122547c6ae99SBarry Smith         }
122647c6ae99SBarry Smith 	ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
122747c6ae99SBarry Smith       }
122847c6ae99SBarry Smith     }
122947c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
123047c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
123147c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
123247c6ae99SBarry Smith   }
123347c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
123447c6ae99SBarry Smith   PetscFunctionReturn(0);
123547c6ae99SBarry Smith }
123647c6ae99SBarry Smith 
123747c6ae99SBarry Smith #undef __FUNCT__
123894013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIBAIJ"
123994013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIBAIJ(DM da,Mat J)
124047c6ae99SBarry Smith {
124147c6ae99SBarry Smith   PetscErrorCode         ierr;
124247c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
124347c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
124447c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
124547c6ae99SBarry Smith   MPI_Comm               comm;
124647c6ae99SBarry Smith   PetscScalar            *values;
12471321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1248aa219208SBarry Smith   DMDAStencilType        st;
124947c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
125047c6ae99SBarry Smith 
125147c6ae99SBarry Smith   PetscFunctionBegin;
125247c6ae99SBarry Smith   /*
125347c6ae99SBarry Smith          nc - number of components per grid point
125447c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
125547c6ae99SBarry Smith 
125647c6ae99SBarry Smith   */
12571321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
125847c6ae99SBarry Smith   col    = 2*s + 1;
125947c6ae99SBarry Smith 
1260aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1261aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
126247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
126347c6ae99SBarry Smith 
126447c6ae99SBarry Smith   ierr  = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
126547c6ae99SBarry Smith 
12661411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
12671411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
126847c6ae99SBarry Smith 
126947c6ae99SBarry Smith   /* determine the matrix preallocation information */
127047c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
127147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
12721321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
12731321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
127447c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
12751321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
12761321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
127747c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
12781321219cSEthan Coon 	kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
12791321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
128047c6ae99SBarry Smith 
128147c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
128247c6ae99SBarry Smith 
128347c6ae99SBarry Smith 	/* Find block columns in block row */
128447c6ae99SBarry Smith 	cnt  = 0;
128547c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
128647c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
128747c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1288aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
128947c6ae99SBarry Smith 		cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
129047c6ae99SBarry Smith 	      }
129147c6ae99SBarry Smith 	    }
129247c6ae99SBarry Smith 	  }
129347c6ae99SBarry Smith 	}
1294784ac674SJed Brown 	ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
129547c6ae99SBarry Smith       }
129647c6ae99SBarry Smith     }
129747c6ae99SBarry Smith   }
129847c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
129947c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
130047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
130147c6ae99SBarry Smith 
1302784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1303784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
130447c6ae99SBarry Smith 
130547c6ae99SBarry Smith   /*
130647c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
130747c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
130847c6ae99SBarry Smith     PETSc ordering.
130947c6ae99SBarry Smith   */
1310fcfd50ebSBarry Smith   if (!da->prealloc_only) {
131147c6ae99SBarry Smith     ierr  = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
131247c6ae99SBarry Smith     ierr  = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
131347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
13141321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
13151321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
131647c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
13171321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
13181321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
131947c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
13201321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
13211321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
132247c6ae99SBarry Smith 
132347c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
132447c6ae99SBarry Smith 
132547c6ae99SBarry Smith 	  cnt  = 0;
132647c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
132747c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
132847c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1329aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
133047c6ae99SBarry Smith                   cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
133147c6ae99SBarry Smith                 }
133247c6ae99SBarry Smith               }
133347c6ae99SBarry Smith             }
133447c6ae99SBarry Smith           }
133547c6ae99SBarry Smith 	  ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
133647c6ae99SBarry Smith 	}
133747c6ae99SBarry Smith       }
133847c6ae99SBarry Smith     }
133947c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
134047c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
134147c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
134247c6ae99SBarry Smith   }
134347c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
134447c6ae99SBarry Smith   PetscFunctionReturn(0);
134547c6ae99SBarry Smith }
134647c6ae99SBarry Smith 
134747c6ae99SBarry Smith #undef __FUNCT__
134847c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular"
134947c6ae99SBarry Smith /*
135047c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
135147c6ae99SBarry Smith   identify in the local ordering with periodic domain.
135247c6ae99SBarry Smith */
135347c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
135447c6ae99SBarry Smith {
135547c6ae99SBarry Smith   PetscErrorCode ierr;
135647c6ae99SBarry Smith   PetscInt       i,n;
135747c6ae99SBarry Smith 
135847c6ae99SBarry Smith   PetscFunctionBegin;
135947c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr);
136047c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr);
136147c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
136247c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
136347c6ae99SBarry Smith   }
136447c6ae99SBarry Smith   *cnt = n;
136547c6ae99SBarry Smith   PetscFunctionReturn(0);
136647c6ae99SBarry Smith }
136747c6ae99SBarry Smith 
136847c6ae99SBarry Smith #undef __FUNCT__
136994013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPISBAIJ"
137094013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPISBAIJ(DM da,Mat J)
137147c6ae99SBarry Smith {
137247c6ae99SBarry Smith   PetscErrorCode         ierr;
137347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
137447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
137547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
137647c6ae99SBarry Smith   MPI_Comm               comm;
137747c6ae99SBarry Smith   PetscScalar            *values;
13781321219cSEthan Coon   DMDABoundaryType       bx,by;
1379aa219208SBarry Smith   DMDAStencilType        st;
138047c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
138147c6ae99SBarry Smith 
138247c6ae99SBarry Smith   PetscFunctionBegin;
138347c6ae99SBarry Smith   /*
138447c6ae99SBarry Smith      nc - number of components per grid point
138547c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
138647c6ae99SBarry Smith   */
13871321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
138847c6ae99SBarry Smith   col = 2*s + 1;
138947c6ae99SBarry Smith 
1390aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1391aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
139247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
139347c6ae99SBarry Smith 
139447c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
139547c6ae99SBarry Smith 
13961411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
13971411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
139847c6ae99SBarry Smith 
139947c6ae99SBarry Smith   /* determine the matrix preallocation information */
140047c6ae99SBarry Smith   ierr = MatPreallocateSymmetricInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
140147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
14021321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14031321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
140447c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
14051321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14061321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
140747c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
140847c6ae99SBarry Smith 
140947c6ae99SBarry Smith       /* Find block columns in block row */
141047c6ae99SBarry Smith       cnt  = 0;
141147c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
141247c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1413aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) {
141447c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
141547c6ae99SBarry Smith           }
141647c6ae99SBarry Smith         }
141747c6ae99SBarry Smith       }
141847c6ae99SBarry Smith       ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
141947c6ae99SBarry Smith       ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
142047c6ae99SBarry Smith     }
142147c6ae99SBarry Smith   }
142247c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
142347c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
142447c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
142547c6ae99SBarry Smith 
1426784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1427784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
142847c6ae99SBarry Smith 
142947c6ae99SBarry Smith   /*
143047c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
143147c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
143247c6ae99SBarry Smith     PETSc ordering.
143347c6ae99SBarry Smith   */
1434fcfd50ebSBarry Smith   if (!da->prealloc_only) {
143547c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
143647c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
143747c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
14381321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14391321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
144047c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
14411321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14421321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
144347c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
144447c6ae99SBarry Smith 
144547c6ae99SBarry Smith         /* Find block columns in block row */
144647c6ae99SBarry Smith         cnt  = 0;
144747c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
144847c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1449aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) {
145047c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
145147c6ae99SBarry Smith             }
145247c6ae99SBarry Smith           }
145347c6ae99SBarry Smith         }
145447c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
145547c6ae99SBarry Smith 	ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
145647c6ae99SBarry Smith       }
145747c6ae99SBarry Smith     }
145847c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
145947c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
146047c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
146147c6ae99SBarry Smith   }
146247c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
146347c6ae99SBarry Smith   PetscFunctionReturn(0);
146447c6ae99SBarry Smith }
146547c6ae99SBarry Smith 
146647c6ae99SBarry Smith #undef __FUNCT__
146794013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPISBAIJ"
146894013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPISBAIJ(DM da,Mat J)
146947c6ae99SBarry Smith {
147047c6ae99SBarry Smith   PetscErrorCode         ierr;
147147c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
147247c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
147347c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
147447c6ae99SBarry Smith   MPI_Comm               comm;
147547c6ae99SBarry Smith   PetscScalar            *values;
14761321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1477aa219208SBarry Smith   DMDAStencilType        st;
147847c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
147947c6ae99SBarry Smith 
148047c6ae99SBarry Smith   PetscFunctionBegin;
148147c6ae99SBarry Smith   /*
148247c6ae99SBarry Smith      nc - number of components per grid point
148347c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
148447c6ae99SBarry Smith   */
14851321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
148647c6ae99SBarry Smith   col = 2*s + 1;
148747c6ae99SBarry Smith 
1488aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1489aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
149047c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
149147c6ae99SBarry Smith 
149247c6ae99SBarry Smith   /* create the matrix */
149347c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
149447c6ae99SBarry Smith 
14951411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
14961411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
149747c6ae99SBarry Smith 
149847c6ae99SBarry Smith   /* determine the matrix preallocation information */
149947c6ae99SBarry Smith   ierr = MatPreallocateSymmetricInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
150047c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
15011321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
15021321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
150347c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
15041321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15051321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
150647c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
15071321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15081321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
150947c6ae99SBarry Smith 
151047c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
151147c6ae99SBarry Smith 
151247c6ae99SBarry Smith 	/* Find block columns in block row */
151347c6ae99SBarry Smith 	cnt  = 0;
151447c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
151547c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
151647c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1517aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
151847c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
151947c6ae99SBarry Smith               }
152047c6ae99SBarry Smith             }
152147c6ae99SBarry Smith           }
152247c6ae99SBarry Smith         }
152347c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
152447c6ae99SBarry Smith         ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
152547c6ae99SBarry Smith       }
152647c6ae99SBarry Smith     }
152747c6ae99SBarry Smith   }
152847c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
152947c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
153047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
153147c6ae99SBarry Smith 
1532784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1533784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
153447c6ae99SBarry Smith 
153547c6ae99SBarry Smith   /*
153647c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
153747c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
153847c6ae99SBarry Smith     PETSc ordering.
153947c6ae99SBarry Smith   */
1540fcfd50ebSBarry Smith   if (!da->prealloc_only) {
154147c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
154247c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));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 	  cnt  = 0;
155647c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
155747c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
155847c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1559aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
156047c6ae99SBarry Smith 		  cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
156147c6ae99SBarry Smith 		}
156247c6ae99SBarry Smith 	      }
156347c6ae99SBarry Smith 	    }
156447c6ae99SBarry Smith 	  }
156547c6ae99SBarry Smith           ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
156647c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
156747c6ae99SBarry Smith 	}
156847c6ae99SBarry Smith       }
156947c6ae99SBarry Smith     }
157047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
157147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157347c6ae99SBarry Smith   }
157447c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
157547c6ae99SBarry Smith   PetscFunctionReturn(0);
157647c6ae99SBarry Smith }
157747c6ae99SBarry Smith 
157847c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
157947c6ae99SBarry Smith 
158047c6ae99SBarry Smith #undef __FUNCT__
158194013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIAIJ_Fill"
158294013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J)
158347c6ae99SBarry Smith {
158447c6ae99SBarry Smith   PetscErrorCode         ierr;
158547c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
158647c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz;
158747c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
158847c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
158947c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
159047c6ae99SBarry Smith   MPI_Comm               comm;
159147c6ae99SBarry Smith   PetscScalar            *values;
15921321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
159347c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1594aa219208SBarry Smith   DMDAStencilType        st;
159547c6ae99SBarry Smith 
159647c6ae99SBarry Smith   PetscFunctionBegin;
159747c6ae99SBarry Smith   /*
159847c6ae99SBarry Smith          nc - number of components per grid point
159947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
160047c6ae99SBarry Smith 
160147c6ae99SBarry Smith   */
16021321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
160347c6ae99SBarry Smith   col    = 2*s + 1;
16041321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){
160547c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
160647c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
160747c6ae99SBarry Smith   }
16081321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){
160947c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
161047c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
161147c6ae99SBarry Smith   }
16121321219cSEthan Coon   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)){
161347c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
161447c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
161547c6ae99SBarry Smith   }
161647c6ae99SBarry Smith 
1617aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1618aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
161947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
162047c6ae99SBarry Smith 
162147c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
16221411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
16231411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
162447c6ae99SBarry Smith 
162547c6ae99SBarry Smith   /* determine the matrix preallocation information */
162647c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
162747c6ae99SBarry Smith 
162847c6ae99SBarry Smith 
162947c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
16301321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16311321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
163247c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
16331321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16341321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
163547c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
16361321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
16371321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
163847c6ae99SBarry Smith 
163947c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
164047c6ae99SBarry Smith 
164147c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
164247c6ae99SBarry Smith 	  cnt  = 0;
164347c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
164447c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
164547c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
164647c6ae99SBarry Smith 		if (ii || jj || kk) {
1647aa219208SBarry Smith 		  if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
164847c6ae99SBarry Smith 		    for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
164947c6ae99SBarry Smith 		      cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
165047c6ae99SBarry Smith 		  }
165147c6ae99SBarry Smith 		} else {
165247c6ae99SBarry Smith 		  if (dfill) {
165347c6ae99SBarry Smith 		    for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
165447c6ae99SBarry Smith 		      cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
165547c6ae99SBarry Smith 		  } else {
165647c6ae99SBarry Smith 		    for (ifill_col=0; ifill_col<nc; ifill_col++)
165747c6ae99SBarry Smith 		      cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
165847c6ae99SBarry Smith 		  }
165947c6ae99SBarry Smith 		}
166047c6ae99SBarry Smith 	      }
166147c6ae99SBarry Smith 	    }
166247c6ae99SBarry Smith 	  }
166347c6ae99SBarry Smith 	  row  = l + nc*(slot);
1664784ac674SJed Brown 	  ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
166547c6ae99SBarry Smith 	}
166647c6ae99SBarry Smith       }
166747c6ae99SBarry Smith     }
166847c6ae99SBarry Smith   }
166947c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
167047c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
167147c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1672784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1673784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
167447c6ae99SBarry Smith 
167547c6ae99SBarry Smith   /*
167647c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
167747c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
167847c6ae99SBarry Smith     PETSc ordering.
167947c6ae99SBarry Smith   */
1680fcfd50ebSBarry Smith   if (!da->prealloc_only) {
168147c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
168247c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
168347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
16841321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16851321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
168647c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
16871321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16881321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
168947c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
16901321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
16911321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
169247c6ae99SBarry Smith 
169347c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
169447c6ae99SBarry Smith 
169547c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
169647c6ae99SBarry Smith 	    cnt  = 0;
169747c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
169847c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
169947c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
170047c6ae99SBarry Smith 		  if (ii || jj || kk) {
1701aa219208SBarry Smith 		    if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
170247c6ae99SBarry Smith 		      for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
170347c6ae99SBarry Smith 			cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
170447c6ae99SBarry Smith 		    }
170547c6ae99SBarry Smith 		  } else {
170647c6ae99SBarry Smith 		    if (dfill) {
170747c6ae99SBarry Smith 		      for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
170847c6ae99SBarry Smith 			cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
170947c6ae99SBarry Smith 		    } else {
171047c6ae99SBarry Smith 		      for (ifill_col=0; ifill_col<nc; ifill_col++)
171147c6ae99SBarry Smith 			cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
171247c6ae99SBarry Smith 		    }
171347c6ae99SBarry Smith 		  }
171447c6ae99SBarry Smith 		}
171547c6ae99SBarry Smith 	      }
171647c6ae99SBarry Smith 	    }
171747c6ae99SBarry Smith 	    row  = l + nc*(slot);
171847c6ae99SBarry Smith 	    ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
171947c6ae99SBarry Smith 	  }
172047c6ae99SBarry Smith 	}
172147c6ae99SBarry Smith       }
172247c6ae99SBarry Smith     }
172347c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
172447c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
172547c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
172647c6ae99SBarry Smith   }
172747c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
172847c6ae99SBarry Smith   PetscFunctionReturn(0);
172947c6ae99SBarry Smith }
1730