xref: /petsc/src/dm/impls/da/fdda.c (revision ce308e1d573627a4a122984c7b5d8134f2a9ebd3)
147c6ae99SBarry Smith 
2b45d2f2cSJed Brown #include <petsc-private/daimpl.h> /*I      "petscdmda.h"     I*/
3c6db04a5SJed Brown #include <petscmat.h>         /*I      "petscmat.h"    I*/
4b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
547c6ae99SBarry Smith 
6e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring *);
7e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring *);
8e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring *);
9e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM,ISColoringType,ISColoring *);
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith /*
1247c6ae99SBarry Smith    For ghost i that may be negative or greater than the upper bound this
1347c6ae99SBarry Smith   maps it into the 0:m-1 range using periodicity
1447c6ae99SBarry Smith */
1547c6ae99SBarry Smith #define SetInRange(i,m) ((i < 0) ? m+i:((i >= m) ? i-m:i))
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith #undef __FUNCT__
18aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills_Private"
19*ce308e1dSBarry Smith static PetscErrorCode DMDASetBlockFills_Private(const PetscInt *dfill,PetscInt w,PetscInt **rfill)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   PetscErrorCode ierr;
2247c6ae99SBarry Smith   PetscInt       i,j,nz,*fill;
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   PetscFunctionBegin;
2547c6ae99SBarry Smith   if (!dfill) PetscFunctionReturn(0);
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   /* count number nonzeros */
2847c6ae99SBarry Smith   nz = 0;
2947c6ae99SBarry Smith   for (i=0; i<w; i++) {
3047c6ae99SBarry Smith     for (j=0; j<w; j++) {
3147c6ae99SBarry Smith       if (dfill[w*i+j]) nz++;
3247c6ae99SBarry Smith     }
3347c6ae99SBarry Smith   }
3447c6ae99SBarry Smith   ierr = PetscMalloc((nz + w + 1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
3547c6ae99SBarry Smith   /* construct modified CSR storage of nonzero structure */
36*ce308e1dSBarry Smith   /*  fill[0 -- w] marks starts of each row of column indices (and end of last row)
37*ce308e1dSBarry Smith    so fill[1] - fill[0] gives number of nonzeros in first row etc */
3847c6ae99SBarry Smith   nz = w + 1;
3947c6ae99SBarry Smith   for (i=0; i<w; i++) {
4047c6ae99SBarry Smith     fill[i] = nz;
4147c6ae99SBarry Smith     for (j=0; j<w; j++) {
4247c6ae99SBarry Smith       if (dfill[w*i+j]) {
4347c6ae99SBarry Smith 	fill[nz] = j;
4447c6ae99SBarry Smith 	nz++;
4547c6ae99SBarry Smith       }
4647c6ae99SBarry Smith     }
4747c6ae99SBarry Smith   }
4847c6ae99SBarry Smith   fill[w] = nz;
4947c6ae99SBarry Smith 
5047c6ae99SBarry Smith   *rfill = fill;
5147c6ae99SBarry Smith   PetscFunctionReturn(0);
5247c6ae99SBarry Smith }
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith #undef __FUNCT__
55aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills"
5647c6ae99SBarry Smith /*@
57aa219208SBarry Smith     DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem
58950540a4SJed Brown     of the matrix returned by DMCreateMatrix().
5947c6ae99SBarry Smith 
60aa219208SBarry Smith     Logically Collective on DMDA
6147c6ae99SBarry Smith 
6247c6ae99SBarry Smith     Input Parameter:
6347c6ae99SBarry Smith +   da - the distributed array
6447c6ae99SBarry Smith .   dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block)
6547c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
6647c6ae99SBarry Smith 
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith     Level: developer
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     Notes: This only makes sense when you are doing multicomponent problems but using the
7147c6ae99SBarry Smith        MPIAIJ matrix format
7247c6ae99SBarry Smith 
7347c6ae99SBarry Smith            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
7447c6ae99SBarry Smith        representing coupling and 0 entries for missing coupling. For example
7547c6ae99SBarry Smith $             dfill[9] = {1, 0, 0,
7647c6ae99SBarry Smith $                         1, 1, 0,
7747c6ae99SBarry Smith $                         0, 1, 1}
7847c6ae99SBarry Smith        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
7947c6ae99SBarry Smith        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
8047c6ae99SBarry Smith        diagonal block).
8147c6ae99SBarry Smith 
82aa219208SBarry Smith      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
8347c6ae99SBarry Smith      can be represented in the dfill, ofill format
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith    Contributed by Glenn Hammond
8647c6ae99SBarry Smith 
878ddb5d8bSBarry Smith .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
8847c6ae99SBarry Smith 
8947c6ae99SBarry Smith @*/
90*ce308e1dSBarry Smith PetscErrorCode  DMDASetBlockFills(DM da,const PetscInt *dfill,const PetscInt *ofill)
9147c6ae99SBarry Smith {
9247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
9347c6ae99SBarry Smith   PetscErrorCode ierr;
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith   PetscFunctionBegin;
96aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr);
97aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr);
9847c6ae99SBarry Smith   PetscFunctionReturn(0);
9947c6ae99SBarry Smith }
10047c6ae99SBarry Smith 
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith #undef __FUNCT__
103e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA"
10419fd82e9SBarry Smith PetscErrorCode  DMCreateColoring_DA(DM da,ISColoringType ctype,MatType mtype,ISColoring *coloring)
10547c6ae99SBarry Smith {
10647c6ae99SBarry Smith   PetscErrorCode   ierr;
10747c6ae99SBarry Smith   PetscInt         dim,m,n,p,nc;
1081321219cSEthan Coon   DMDABoundaryType bx,by,bz;
10947c6ae99SBarry Smith   MPI_Comm         comm;
11047c6ae99SBarry Smith   PetscMPIInt      size;
11147c6ae99SBarry Smith   PetscBool        isBAIJ;
11247c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
11347c6ae99SBarry Smith 
11447c6ae99SBarry Smith   PetscFunctionBegin;
11547c6ae99SBarry Smith   /*
11647c6ae99SBarry Smith                                   m
11747c6ae99SBarry Smith           ------------------------------------------------------
11847c6ae99SBarry Smith          |                                                     |
11947c6ae99SBarry Smith          |                                                     |
12047c6ae99SBarry Smith          |               ----------------------                |
12147c6ae99SBarry Smith          |               |                    |                |
12247c6ae99SBarry Smith       n  |           yn  |                    |                |
12347c6ae99SBarry Smith          |               |                    |                |
12447c6ae99SBarry Smith          |               .---------------------                |
12547c6ae99SBarry Smith          |             (xs,ys)     xn                          |
12647c6ae99SBarry Smith          |            .                                        |
12747c6ae99SBarry Smith          |         (gxs,gys)                                   |
12847c6ae99SBarry Smith          |                                                     |
12947c6ae99SBarry Smith           -----------------------------------------------------
13047c6ae99SBarry Smith   */
13147c6ae99SBarry Smith 
13247c6ae99SBarry Smith   /*
13347c6ae99SBarry Smith          nc - number of components per grid point
13447c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
13547c6ae99SBarry Smith 
13647c6ae99SBarry Smith   */
1371321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&bx,&by,&bz,0);CHKERRQ(ierr);
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
14047c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
14147c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED){
14247c6ae99SBarry Smith     if (size == 1) {
14347c6ae99SBarry Smith       ctype = IS_COLORING_GLOBAL;
14447c6ae99SBarry Smith     } else if (dim > 1){
1451321219cSEthan Coon       if ((m==1 && bx == DMDA_BOUNDARY_PERIODIC) || (n==1 && by == DMDA_BOUNDARY_PERIODIC) || (p==1 && bz == DMDA_BOUNDARY_PERIODIC)){
14647c6ae99SBarry 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");
14747c6ae99SBarry Smith       }
14847c6ae99SBarry Smith     }
14947c6ae99SBarry Smith   }
15047c6ae99SBarry Smith 
151aa219208SBarry Smith   /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ
15247c6ae99SBarry Smith      matrices is for the blocks, not the individual matrix elements  */
1534833614aSSatish Balay   ierr = PetscStrcmp(mtype,MATBAIJ,&isBAIJ);CHKERRQ(ierr);
1544833614aSSatish Balay   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);}
15547c6ae99SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);}
15647c6ae99SBarry Smith   if (isBAIJ) {
15747c6ae99SBarry Smith     dd->w = 1;
15847c6ae99SBarry Smith     dd->xs = dd->xs/nc;
15947c6ae99SBarry Smith     dd->xe = dd->xe/nc;
16047c6ae99SBarry Smith     dd->Xs = dd->Xs/nc;
16147c6ae99SBarry Smith     dd->Xe = dd->Xe/nc;
16247c6ae99SBarry Smith   }
16347c6ae99SBarry Smith 
16447c6ae99SBarry Smith   /*
165aa219208SBarry Smith      We do not provide a getcoloring function in the DMDA operations because
166aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
16747c6ae99SBarry Smith    more low-level then matrices.
16847c6ae99SBarry Smith   */
16947c6ae99SBarry Smith   if (dim == 1) {
170e727c939SJed Brown     ierr = DMCreateColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
17147c6ae99SBarry Smith   } else if (dim == 2) {
172e727c939SJed Brown     ierr =  DMCreateColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
17347c6ae99SBarry Smith   } else if (dim == 3) {
174e727c939SJed Brown     ierr =  DMCreateColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
17571cd77b2SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim);
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__
189e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_MPIAIJ"
190e727c939SJed Brown PetscErrorCode DMCreateColoring_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) {
215e727c939SJed Brown     ierr = DMCreateColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
21647c6ae99SBarry Smith   } else {
21747c6ae99SBarry Smith 
21866a15934SBarry Smith     if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\
21947c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", m, col);
22066a15934SBarry Smith     if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\
22147c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", n, col);
22247c6ae99SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
22347c6ae99SBarry Smith       if (!dd->localcoloring) {
22447c6ae99SBarry Smith 	ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
22547c6ae99SBarry Smith 	ii = 0;
22647c6ae99SBarry Smith 	for (j=ys; j<ys+ny; j++) {
22747c6ae99SBarry Smith 	  for (i=xs; i<xs+nx; i++) {
22847c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
22947c6ae99SBarry Smith 	      colors[ii++] = k + nc*((i % col) + col*(j % col));
23047c6ae99SBarry Smith 	    }
23147c6ae99SBarry Smith 	  }
23247c6ae99SBarry Smith 	}
23347c6ae99SBarry Smith         ncolors = nc + nc*(col-1 + col*(col-1));
23447c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
23547c6ae99SBarry Smith       }
23647c6ae99SBarry Smith       *coloring = dd->localcoloring;
23747c6ae99SBarry Smith     } else if (ctype == IS_COLORING_GHOSTED) {
23847c6ae99SBarry Smith       if (!dd->ghostedcoloring) {
23947c6ae99SBarry Smith 	ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
24047c6ae99SBarry Smith 	ii = 0;
24147c6ae99SBarry Smith 	for (j=gys; j<gys+gny; j++) {
24247c6ae99SBarry Smith 	  for (i=gxs; i<gxs+gnx; i++) {
24347c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
24447c6ae99SBarry Smith 	      /* the complicated stuff is to handle periodic boundaries */
24547c6ae99SBarry Smith 	      colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
24647c6ae99SBarry Smith 	    }
24747c6ae99SBarry Smith 	  }
24847c6ae99SBarry Smith 	}
24947c6ae99SBarry Smith         ncolors = nc + nc*(col - 1 + col*(col-1));
25047c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
25147c6ae99SBarry Smith         /* PetscIntView(ncolors,(PetscInt *)colors,0); */
25247c6ae99SBarry Smith 
25347c6ae99SBarry Smith 	ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
25447c6ae99SBarry Smith       }
25547c6ae99SBarry Smith       *coloring = dd->ghostedcoloring;
25647c6ae99SBarry Smith     } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
25747c6ae99SBarry Smith   }
25847c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
25947c6ae99SBarry Smith   PetscFunctionReturn(0);
26047c6ae99SBarry Smith }
26147c6ae99SBarry Smith 
26247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
26347c6ae99SBarry Smith 
26447c6ae99SBarry Smith #undef __FUNCT__
265e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_3d_MPIAIJ"
266e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
26747c6ae99SBarry Smith {
26847c6ae99SBarry Smith   PetscErrorCode    ierr;
26947c6ae99SBarry 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;
27047c6ae99SBarry Smith   PetscInt          ncolors;
27147c6ae99SBarry Smith   MPI_Comm          comm;
2721321219cSEthan Coon   DMDABoundaryType  bx,by,bz;
273aa219208SBarry Smith   DMDAStencilType   st;
27447c6ae99SBarry Smith   ISColoringValue   *colors;
27547c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
27647c6ae99SBarry Smith 
27747c6ae99SBarry Smith   PetscFunctionBegin;
27847c6ae99SBarry Smith   /*
27947c6ae99SBarry Smith          nc - number of components per grid point
28047c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
28147c6ae99SBarry Smith 
28247c6ae99SBarry Smith   */
2831321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
28447c6ae99SBarry Smith   col    = 2*s + 1;
28566a15934SBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
28647c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
28766a15934SBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
28847c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
28966a15934SBarry Smith   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
29047c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
29147c6ae99SBarry Smith 
292aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
293aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
29447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
29547c6ae99SBarry Smith 
29647c6ae99SBarry Smith   /* create the coloring */
29747c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
29847c6ae99SBarry Smith     if (!dd->localcoloring) {
29947c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
30047c6ae99SBarry Smith       ii = 0;
30147c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
30247c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
30347c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
30447c6ae99SBarry Smith             for (l=0; l<nc; l++) {
30547c6ae99SBarry Smith               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
30647c6ae99SBarry Smith             }
30747c6ae99SBarry Smith           }
30847c6ae99SBarry Smith         }
30947c6ae99SBarry Smith       }
31047c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
31147c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr);
31247c6ae99SBarry Smith     }
31347c6ae99SBarry Smith     *coloring = dd->localcoloring;
31447c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
31547c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
31647c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
31747c6ae99SBarry Smith       ii = 0;
31847c6ae99SBarry Smith       for (k=gzs; k<gzs+gnz; k++) {
31947c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
32047c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
32147c6ae99SBarry Smith             for (l=0; l<nc; l++) {
32247c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
32347c6ae99SBarry Smith               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
32447c6ae99SBarry Smith             }
32547c6ae99SBarry Smith           }
32647c6ae99SBarry Smith         }
32747c6ae99SBarry Smith       }
32847c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
32947c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
33047c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
33147c6ae99SBarry Smith     }
33247c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
33347c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
33447c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
33547c6ae99SBarry Smith   PetscFunctionReturn(0);
33647c6ae99SBarry Smith }
33747c6ae99SBarry Smith 
33847c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
33947c6ae99SBarry Smith 
34047c6ae99SBarry Smith #undef __FUNCT__
341e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_1d_MPIAIJ"
342e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
34347c6ae99SBarry Smith {
34447c6ae99SBarry Smith   PetscErrorCode    ierr;
34547c6ae99SBarry Smith   PetscInt          xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
34647c6ae99SBarry Smith   PetscInt          ncolors;
34747c6ae99SBarry Smith   MPI_Comm          comm;
3481321219cSEthan Coon   DMDABoundaryType  bx;
34947c6ae99SBarry Smith   ISColoringValue   *colors;
35047c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith   PetscFunctionBegin;
35347c6ae99SBarry Smith   /*
35447c6ae99SBarry Smith          nc - number of components per grid point
35547c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
35647c6ae99SBarry Smith 
35747c6ae99SBarry Smith   */
3581321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
35947c6ae99SBarry Smith   col    = 2*s + 1;
36047c6ae99SBarry Smith 
36166a15934SBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points %d is divisible\n\
36231e6f798SBarry Smith                                                           by 2*stencil_width + 1 %d\n",(int)m,(int)col);
36347c6ae99SBarry Smith 
364aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
365aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
36647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
36747c6ae99SBarry Smith 
36847c6ae99SBarry Smith   /* create the coloring */
36947c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
37047c6ae99SBarry Smith     if (!dd->localcoloring) {
37147c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
37247c6ae99SBarry Smith       i1 = 0;
37347c6ae99SBarry Smith       for (i=xs; i<xs+nx; i++) {
37447c6ae99SBarry Smith         for (l=0; l<nc; l++) {
37547c6ae99SBarry Smith           colors[i1++] = l + nc*(i % col);
37647c6ae99SBarry Smith         }
37747c6ae99SBarry Smith       }
37847c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
37947c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr);
38047c6ae99SBarry Smith     }
38147c6ae99SBarry Smith     *coloring = dd->localcoloring;
38247c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
38347c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
38447c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
38547c6ae99SBarry Smith       i1 = 0;
38647c6ae99SBarry Smith       for (i=gxs; i<gxs+gnx; i++) {
38747c6ae99SBarry Smith         for (l=0; l<nc; l++) {
38847c6ae99SBarry Smith           /* the complicated stuff is to handle periodic boundaries */
38947c6ae99SBarry Smith           colors[i1++] = l + nc*(SetInRange(i,m) % col);
39047c6ae99SBarry Smith         }
39147c6ae99SBarry Smith       }
39247c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
39347c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
39447c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
39547c6ae99SBarry Smith     }
39647c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
39747c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
39847c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
39947c6ae99SBarry Smith   PetscFunctionReturn(0);
40047c6ae99SBarry Smith }
40147c6ae99SBarry Smith 
40247c6ae99SBarry Smith #undef __FUNCT__
403e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_5pt_MPIAIJ"
404e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
40547c6ae99SBarry Smith {
40647c6ae99SBarry Smith   PetscErrorCode    ierr;
40747c6ae99SBarry Smith   PetscInt          xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
40847c6ae99SBarry Smith   PetscInt          ncolors;
40947c6ae99SBarry Smith   MPI_Comm          comm;
4101321219cSEthan Coon   DMDABoundaryType  bx,by;
41147c6ae99SBarry Smith   ISColoringValue   *colors;
41247c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
41347c6ae99SBarry Smith 
41447c6ae99SBarry Smith   PetscFunctionBegin;
41547c6ae99SBarry Smith   /*
41647c6ae99SBarry Smith          nc - number of components per grid point
41747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
41847c6ae99SBarry Smith 
41947c6ae99SBarry Smith   */
4201321219cSEthan Coon   ierr   = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);CHKERRQ(ierr);
421aa219208SBarry Smith   ierr   = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
422aa219208SBarry Smith   ierr   = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
42347c6ae99SBarry Smith   ierr   = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
42447c6ae99SBarry Smith 
42566a15934SBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % 5)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible by 5\n");
42666a15934SBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % 5)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible by 5\n");
42747c6ae99SBarry Smith 
42847c6ae99SBarry Smith   /* create the coloring */
42947c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
43047c6ae99SBarry Smith     if (!dd->localcoloring) {
43147c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
43247c6ae99SBarry Smith       ii = 0;
43347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
43447c6ae99SBarry Smith 	for (i=xs; i<xs+nx; i++) {
43547c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
43647c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*j+i) % 5);
43747c6ae99SBarry Smith 	  }
43847c6ae99SBarry Smith 	}
43947c6ae99SBarry Smith       }
44047c6ae99SBarry Smith       ncolors = 5*nc;
44147c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
44247c6ae99SBarry Smith     }
44347c6ae99SBarry Smith     *coloring = dd->localcoloring;
44447c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
44547c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
44647c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
44747c6ae99SBarry Smith       ii = 0;
44847c6ae99SBarry Smith       for (j=gys; j<gys+gny; j++) {
44947c6ae99SBarry Smith 	for (i=gxs; i<gxs+gnx; i++) {
45047c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
45147c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
45247c6ae99SBarry Smith 	  }
45347c6ae99SBarry Smith 	}
45447c6ae99SBarry Smith       }
45547c6ae99SBarry Smith       ncolors = 5*nc;
45647c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
45747c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
45847c6ae99SBarry Smith     }
45947c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
46047c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
46147c6ae99SBarry Smith   PetscFunctionReturn(0);
46247c6ae99SBarry Smith }
46347c6ae99SBarry Smith 
46447c6ae99SBarry Smith /* =========================================================================== */
465950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM,Mat);
466*ce308e1dSBarry Smith extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM,Mat);
467950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM,Mat);
468950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM,Mat);
469950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM,Mat);
470950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM,Mat);
471950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM,Mat);
472950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM,Mat);
473950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM,Mat);
474950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM,Mat);
47547c6ae99SBarry Smith 
47647c6ae99SBarry Smith #undef __FUNCT__
477c688c046SMatthew G Knepley #define __FUNCT__ "MatSetupDM"
4788bbdbebaSMatthew G Knepley /*@C
479c688c046SMatthew G Knepley    MatSetupDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix
48047c6ae99SBarry Smith 
48147c6ae99SBarry Smith    Logically Collective on Mat
48247c6ae99SBarry Smith 
48347c6ae99SBarry Smith    Input Parameters:
48447c6ae99SBarry Smith +  mat - the matrix
48547c6ae99SBarry Smith -  da - the da
48647c6ae99SBarry Smith 
48747c6ae99SBarry Smith    Level: intermediate
48847c6ae99SBarry Smith 
48947c6ae99SBarry Smith @*/
490c688c046SMatthew G Knepley PetscErrorCode MatSetupDM(Mat mat,DM da)
49147c6ae99SBarry Smith {
49247c6ae99SBarry Smith   PetscErrorCode ierr;
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith   PetscFunctionBegin;
49547c6ae99SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
49647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
497c688c046SMatthew G Knepley   ierr = PetscTryMethod(mat,"MatSetupDM_C",(Mat,DM),(mat,da));CHKERRQ(ierr);
49847c6ae99SBarry Smith   PetscFunctionReturn(0);
49947c6ae99SBarry Smith }
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith EXTERN_C_BEGIN
50247c6ae99SBarry Smith #undef __FUNCT__
50347c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA"
5047087cfbeSBarry Smith PetscErrorCode  MatView_MPI_DA(Mat A,PetscViewer viewer)
50547c6ae99SBarry Smith {
5069a42bb27SBarry Smith   DM             da;
50747c6ae99SBarry Smith   PetscErrorCode ierr;
50847c6ae99SBarry Smith   const char     *prefix;
50947c6ae99SBarry Smith   Mat            Anatural;
51047c6ae99SBarry Smith   AO             ao;
51147c6ae99SBarry Smith   PetscInt       rstart,rend,*petsc,i;
51247c6ae99SBarry Smith   IS             is;
51347c6ae99SBarry Smith   MPI_Comm       comm;
51474388724SJed Brown   PetscViewerFormat format;
51547c6ae99SBarry Smith 
51647c6ae99SBarry Smith   PetscFunctionBegin;
51774388724SJed Brown   /* Check whether we are just printing info, in which case MatView() already viewed everything we wanted to view */
51874388724SJed Brown   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
51974388724SJed Brown   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
52074388724SJed Brown 
52147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
522c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
523aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
52447c6ae99SBarry Smith 
525aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
52647c6ae99SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
52747c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr);
52847c6ae99SBarry Smith   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
52947c6ae99SBarry Smith   ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr);
53047c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
53147c6ae99SBarry Smith 
53247c6ae99SBarry Smith   /* call viewer on natural ordering */
53347c6ae99SBarry Smith   ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr);
534fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
53547c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr);
53647c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr);
53747c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr);
53847c6ae99SBarry Smith   ierr = MatView(Anatural,viewer);CHKERRQ(ierr);
539fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
54047c6ae99SBarry Smith   PetscFunctionReturn(0);
54147c6ae99SBarry Smith }
54247c6ae99SBarry Smith EXTERN_C_END
54347c6ae99SBarry Smith 
54447c6ae99SBarry Smith EXTERN_C_BEGIN
54547c6ae99SBarry Smith #undef __FUNCT__
54647c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA"
5477087cfbeSBarry Smith PetscErrorCode  MatLoad_MPI_DA(Mat A,PetscViewer viewer)
54847c6ae99SBarry Smith {
5499a42bb27SBarry Smith   DM             da;
55047c6ae99SBarry Smith   PetscErrorCode ierr;
55147c6ae99SBarry Smith   Mat            Anatural,Aapp;
55247c6ae99SBarry Smith   AO             ao;
55347c6ae99SBarry Smith   PetscInt       rstart,rend,*app,i;
55447c6ae99SBarry Smith   IS             is;
55547c6ae99SBarry Smith   MPI_Comm       comm;
55647c6ae99SBarry Smith 
55747c6ae99SBarry Smith   PetscFunctionBegin;
55847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
559c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
560aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith   /* Load the matrix in natural ordering */
56347c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,&Anatural);CHKERRQ(ierr);
56447c6ae99SBarry Smith   ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr);
56547c6ae99SBarry Smith   ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
56647c6ae99SBarry Smith   ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr);
56747c6ae99SBarry Smith 
56847c6ae99SBarry Smith   /* Map natural ordering to application ordering and create IS */
569aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
57047c6ae99SBarry Smith   ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr);
57147c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr);
57247c6ae99SBarry Smith   for (i=rstart; i<rend; i++) app[i-rstart] = i;
57347c6ae99SBarry Smith   ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr);
57447c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
57547c6ae99SBarry Smith 
57647c6ae99SBarry Smith   /* Do permutation and replace header */
57747c6ae99SBarry Smith   ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr);
57847c6ae99SBarry Smith   ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr);
579fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
580fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
58147c6ae99SBarry Smith   PetscFunctionReturn(0);
58247c6ae99SBarry Smith }
58347c6ae99SBarry Smith EXTERN_C_END
58447c6ae99SBarry Smith 
58547c6ae99SBarry Smith #undef __FUNCT__
586950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA"
58719fd82e9SBarry Smith PetscErrorCode DMCreateMatrix_DA(DM da, MatType mtype,Mat *J)
58847c6ae99SBarry Smith {
58947c6ae99SBarry Smith   PetscErrorCode ierr;
59047c6ae99SBarry Smith   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
59147c6ae99SBarry Smith   Mat            A;
59247c6ae99SBarry Smith   MPI_Comm       comm;
59319fd82e9SBarry Smith   MatType        Atype;
59437d0c07bSMatthew G Knepley   PetscSection   section, sectionGlobal;
59547c6ae99SBarry Smith   void           (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL;
59647c6ae99SBarry Smith   MatType        ttype[256];
59747c6ae99SBarry Smith   PetscBool      flg;
59847c6ae99SBarry Smith   PetscMPIInt    size;
59947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
60047c6ae99SBarry Smith 
60147c6ae99SBarry Smith   PetscFunctionBegin;
60247c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
60347c6ae99SBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
60447c6ae99SBarry Smith #endif
6055da5aae0SJed Brown   if (!mtype) mtype = MATAIJ;
60647c6ae99SBarry Smith   ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr);
607aa219208SBarry Smith   ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DMDA options","Mat");CHKERRQ(ierr);
608dd85299cSBarry Smith   ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr);
60947c6ae99SBarry Smith   ierr = PetscOptionsEnd();
61047c6ae99SBarry Smith 
61137d0c07bSMatthew G Knepley   ierr = DMGetDefaultSection(da, &section);CHKERRQ(ierr);
61237d0c07bSMatthew G Knepley   if (section) {
61337d0c07bSMatthew G Knepley     PetscInt  bs = -1;
61437d0c07bSMatthew G Knepley     PetscInt  localSize;
61537d0c07bSMatthew G Knepley     PetscBool isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isSymmetric;
61637d0c07bSMatthew G Knepley 
61737d0c07bSMatthew G Knepley     ierr = DMGetDefaultGlobalSection(da, &sectionGlobal);CHKERRQ(ierr);
61837d0c07bSMatthew G Knepley     ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize);CHKERRQ(ierr);
61937d0c07bSMatthew G Knepley     ierr = MatCreate(((PetscObject) da)->comm, J);CHKERRQ(ierr);
62037d0c07bSMatthew G Knepley     ierr = MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
62137d0c07bSMatthew G Knepley     ierr = MatSetType(*J, mtype);CHKERRQ(ierr);
62237d0c07bSMatthew G Knepley     ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
62337d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSHELL, &isShell);CHKERRQ(ierr);
62437d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATBAIJ, &isBlock);CHKERRQ(ierr);
62537d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSEQBAIJ, &isSeqBlock);CHKERRQ(ierr);
62637d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATMPIBAIJ, &isMPIBlock);CHKERRQ(ierr);
62737d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr);
62837d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr);
62937d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr);
63037d0c07bSMatthew G Knepley     /* Check for symmetric storage */
63137d0c07bSMatthew G Knepley     isSymmetric = (PetscBool) (isSymBlock || isSymSeqBlock || isSymMPIBlock);
63237d0c07bSMatthew G Knepley     if (isSymmetric) {
63337d0c07bSMatthew G Knepley       ierr = MatSetOption(*J, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);
63437d0c07bSMatthew G Knepley     }
63537d0c07bSMatthew G Knepley     if (!isShell) {
6364020307fSSatish Balay       /* PetscBool fillMatrix = (PetscBool) !da->prealloc_only; */
63737d0c07bSMatthew G Knepley       PetscInt *dnz, *onz, *dnzu, *onzu, bsLocal;
63837d0c07bSMatthew G Knepley 
63937d0c07bSMatthew G Knepley       if (bs < 0) {
64037d0c07bSMatthew G Knepley         if (isBlock || isSeqBlock || isMPIBlock || isSymBlock || isSymSeqBlock || isSymMPIBlock) {
64137d0c07bSMatthew G Knepley           PetscInt pStart, pEnd, p, dof;
64237d0c07bSMatthew G Knepley 
64337d0c07bSMatthew G Knepley           ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
64437d0c07bSMatthew G Knepley           for (p = pStart; p < pEnd; ++p) {
64537d0c07bSMatthew G Knepley             ierr = PetscSectionGetDof(sectionGlobal, p, &dof);CHKERRQ(ierr);
64637d0c07bSMatthew G Knepley             if (dof) {
64737d0c07bSMatthew G Knepley               bs = dof;
64837d0c07bSMatthew G Knepley               break;
64937d0c07bSMatthew G Knepley             }
65037d0c07bSMatthew G Knepley           }
65137d0c07bSMatthew G Knepley         } else {
65237d0c07bSMatthew G Knepley           bs = 1;
65337d0c07bSMatthew G Knepley         }
65437d0c07bSMatthew G Knepley         /* Must have same blocksize on all procs (some might have no points) */
65537d0c07bSMatthew G Knepley         bsLocal = bs;
65637d0c07bSMatthew G Knepley         ierr = MPI_Allreduce(&bsLocal, &bs, 1, MPIU_INT, MPI_MAX, ((PetscObject) da)->comm);CHKERRQ(ierr);
65737d0c07bSMatthew G Knepley       }
65837d0c07bSMatthew G Knepley       ierr = PetscMalloc4(localSize/bs, PetscInt, &dnz, localSize/bs, PetscInt, &onz, localSize/bs, PetscInt, &dnzu, localSize/bs, PetscInt, &onzu);CHKERRQ(ierr);
65937d0c07bSMatthew G Knepley       ierr = PetscMemzero(dnz,  localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
66037d0c07bSMatthew G Knepley       ierr = PetscMemzero(onz,  localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
66137d0c07bSMatthew G Knepley       ierr = PetscMemzero(dnzu, localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
66237d0c07bSMatthew G Knepley       ierr = PetscMemzero(onzu, localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
66337d0c07bSMatthew G Knepley       /* ierr = DMComplexPreallocateOperator(dm, bs, section, sectionGlobal, dnz, onz, dnzu, onzu, *J, fillMatrix);CHKERRQ(ierr); */
66437d0c07bSMatthew G Knepley       ierr = PetscFree4(dnz, onz, dnzu, onzu);CHKERRQ(ierr);
66537d0c07bSMatthew G Knepley     }
66637d0c07bSMatthew G Knepley   }
66747c6ae99SBarry Smith   /*
66847c6ae99SBarry Smith                                   m
66947c6ae99SBarry Smith           ------------------------------------------------------
67047c6ae99SBarry Smith          |                                                     |
67147c6ae99SBarry Smith          |                                                     |
67247c6ae99SBarry Smith          |               ----------------------                |
67347c6ae99SBarry Smith          |               |                    |                |
67447c6ae99SBarry Smith       n  |           ny  |                    |                |
67547c6ae99SBarry Smith          |               |                    |                |
67647c6ae99SBarry Smith          |               .---------------------                |
67747c6ae99SBarry Smith          |             (xs,ys)     nx                          |
67847c6ae99SBarry Smith          |            .                                        |
67947c6ae99SBarry Smith          |         (gxs,gys)                                   |
68047c6ae99SBarry Smith          |                                                     |
68147c6ae99SBarry Smith           -----------------------------------------------------
68247c6ae99SBarry Smith   */
68347c6ae99SBarry Smith 
68447c6ae99SBarry Smith   /*
68547c6ae99SBarry Smith          nc - number of components per grid point
68647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
68747c6ae99SBarry Smith 
68847c6ae99SBarry Smith   */
689e30e807fSPeter Brune   M = dd->M;
690e30e807fSPeter Brune   N = dd->N;
691e30e807fSPeter Brune   P = dd->P;
692e30e807fSPeter Brune   dim = dd->dim;
693e30e807fSPeter Brune   dof = dd->w;
694e30e807fSPeter Brune   /* ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); */
695aa219208SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
69647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
69747c6ae99SBarry Smith   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
69847c6ae99SBarry Smith   ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr);
69919fd82e9SBarry Smith   ierr = MatSetType(A,(MatType)ttype);CHKERRQ(ierr);
70095ee5b0eSBarry Smith   ierr = MatSetDM(A,da);CHKERRQ(ierr);
70147c6ae99SBarry Smith   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
70247c6ae99SBarry Smith   ierr = MatGetType(A,&Atype);CHKERRQ(ierr);
70347c6ae99SBarry Smith   /*
704aa219208SBarry Smith      We do not provide a getmatrix function in the DMDA operations because
705aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
70647c6ae99SBarry Smith    more low-level than matrices. This is kind of cheating but, cause sometimes
707aa219208SBarry Smith    we think of DMDA has higher level than matrices.
70847c6ae99SBarry Smith 
70947c6ae99SBarry Smith      We could switch based on Atype (or mtype), but we do not since the
71047c6ae99SBarry Smith    specialized setting routines depend only the particular preallocation
71147c6ae99SBarry Smith    details of the matrix, not the type itself.
71247c6ae99SBarry Smith   */
71347c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
71447c6ae99SBarry Smith   if (!aij) {
71547c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
71647c6ae99SBarry Smith   }
71747c6ae99SBarry Smith   if (!aij) {
71847c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
71947c6ae99SBarry Smith     if (!baij) {
72047c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
72147c6ae99SBarry Smith     }
72247c6ae99SBarry Smith     if (!baij){
72347c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
72447c6ae99SBarry Smith       if (!sbaij) {
72547c6ae99SBarry Smith         ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
72647c6ae99SBarry Smith       }
72747c6ae99SBarry Smith     }
72847c6ae99SBarry Smith   }
72947c6ae99SBarry Smith   if (aij) {
73047c6ae99SBarry Smith     if (dim == 1) {
731*ce308e1dSBarry Smith       if (dd->ofill) {
732*ce308e1dSBarry Smith         ierr = DMCreateMatrix_DA_1d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
733*ce308e1dSBarry Smith       } else {
734950540a4SJed Brown         ierr = DMCreateMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr);
735*ce308e1dSBarry Smith       }
73647c6ae99SBarry Smith     } else if (dim == 2) {
73747c6ae99SBarry Smith       if (dd->ofill) {
738950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
73947c6ae99SBarry Smith       } else {
740950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr);
74147c6ae99SBarry Smith       }
74247c6ae99SBarry Smith     } else if (dim == 3) {
74347c6ae99SBarry Smith       if (dd->ofill) {
744950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
74547c6ae99SBarry Smith       } else {
746950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr);
74747c6ae99SBarry Smith       }
74847c6ae99SBarry Smith     }
74947c6ae99SBarry Smith   } else if (baij) {
75047c6ae99SBarry Smith     if (dim == 2) {
751950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr);
75247c6ae99SBarry Smith     } else if (dim == 3) {
753950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr);
75466a15934SBarry Smith     } else  SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension! Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim);
75547c6ae99SBarry Smith   } else if (sbaij) {
75647c6ae99SBarry Smith     if (dim == 2) {
757950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr);
75847c6ae99SBarry Smith     } else if (dim == 3) {
759950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr);
76066a15934SBarry Smith     } else SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension! Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim);
761869776cdSLisandro Dalcin   } else {
762869776cdSLisandro Dalcin     ISLocalToGlobalMapping ltog,ltogb;
763869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
764869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
7652949035bSJed Brown     ierr = MatSetUp(A);CHKERRQ(ierr);
766869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMapping(A,ltog,ltog);CHKERRQ(ierr);
767869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMappingBlock(A,ltogb,ltogb);CHKERRQ(ierr);
76847c6ae99SBarry Smith   }
769aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
77047c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
771c688c046SMatthew G Knepley   ierr = MatSetDM(A,da);CHKERRQ(ierr);
77247c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
77347c6ae99SBarry Smith   if (size > 1) {
77447c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
77547c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);CHKERRQ(ierr);
77647c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);CHKERRQ(ierr);
77747c6ae99SBarry Smith   }
77847c6ae99SBarry Smith   *J = A;
77947c6ae99SBarry Smith   PetscFunctionReturn(0);
78047c6ae99SBarry Smith }
78147c6ae99SBarry Smith 
78247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
78347c6ae99SBarry Smith #undef __FUNCT__
784950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ"
785950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM da,Mat J)
78647c6ae99SBarry Smith {
78747c6ae99SBarry Smith   PetscErrorCode         ierr;
78847c6ae99SBarry 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;
78947c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
79047c6ae99SBarry Smith   MPI_Comm               comm;
79147c6ae99SBarry Smith   PetscScalar            *values;
7921321219cSEthan Coon   DMDABoundaryType       bx,by;
79347c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
794aa219208SBarry Smith   DMDAStencilType        st;
79547c6ae99SBarry Smith 
79647c6ae99SBarry Smith   PetscFunctionBegin;
79747c6ae99SBarry Smith   /*
79847c6ae99SBarry Smith          nc - number of components per grid point
79947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
80047c6ae99SBarry Smith 
80147c6ae99SBarry Smith   */
8021321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
80347c6ae99SBarry Smith   col = 2*s + 1;
804aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
805aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
80647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
80747c6ae99SBarry Smith 
80847c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
8091411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
8101411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
81147c6ae99SBarry Smith 
81247c6ae99SBarry Smith   /* determine the matrix preallocation information */
81347c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
81447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
81547c6ae99SBarry Smith 
8161321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8171321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
81847c6ae99SBarry Smith 
81947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
82047c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
82147c6ae99SBarry Smith 
8221321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8231321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
82447c6ae99SBarry Smith 
82547c6ae99SBarry Smith       cnt  = 0;
82647c6ae99SBarry Smith       for (k=0; k<nc; k++) {
82747c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
82847c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
829aa219208SBarry Smith 	    if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
83047c6ae99SBarry Smith 	      cols[cnt++]  = k + nc*(slot + gnx*l + p);
83147c6ae99SBarry Smith 	    }
83247c6ae99SBarry Smith 	  }
83347c6ae99SBarry Smith 	}
83447c6ae99SBarry Smith 	rows[k] = k + nc*(slot);
83547c6ae99SBarry Smith       }
836784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
83747c6ae99SBarry Smith     }
83847c6ae99SBarry Smith   }
839f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
84047c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
84147c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
84247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
84347c6ae99SBarry Smith 
844784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
845784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
84647c6ae99SBarry Smith 
84747c6ae99SBarry Smith   /*
84847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
84947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
85047c6ae99SBarry Smith     PETSc ordering.
85147c6ae99SBarry Smith   */
852fcfd50ebSBarry Smith   if (!da->prealloc_only) {
85347c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
85447c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
85547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
85647c6ae99SBarry Smith 
8571321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8581321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
86147c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
86247c6ae99SBarry Smith 
8631321219cSEthan Coon 	lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8641321219cSEthan Coon 	lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
86547c6ae99SBarry Smith 
86647c6ae99SBarry Smith 	cnt  = 0;
86747c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
86847c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
86947c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
870aa219208SBarry Smith 	      if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
87147c6ae99SBarry Smith 		cols[cnt++]  = k + nc*(slot + gnx*l + p);
87247c6ae99SBarry Smith 	      }
87347c6ae99SBarry Smith 	    }
87447c6ae99SBarry Smith 	  }
87547c6ae99SBarry Smith 	  rows[k]      = k + nc*(slot);
87647c6ae99SBarry Smith 	}
87747c6ae99SBarry Smith 	ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
87847c6ae99SBarry Smith       }
87947c6ae99SBarry Smith     }
88047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
88147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
88247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
88347c6ae99SBarry Smith   }
88447c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
88547c6ae99SBarry Smith   PetscFunctionReturn(0);
88647c6ae99SBarry Smith }
88747c6ae99SBarry Smith 
88847c6ae99SBarry Smith #undef __FUNCT__
889950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ_Fill"
890950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J)
89147c6ae99SBarry Smith {
89247c6ae99SBarry Smith   PetscErrorCode         ierr;
89347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
89447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p;
89547c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
89647c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
89747c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
89847c6ae99SBarry Smith   MPI_Comm               comm;
89947c6ae99SBarry Smith   PetscScalar            *values;
9001321219cSEthan Coon   DMDABoundaryType       bx,by;
90147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
902aa219208SBarry Smith   DMDAStencilType        st;
90347c6ae99SBarry Smith 
90447c6ae99SBarry Smith   PetscFunctionBegin;
90547c6ae99SBarry Smith   /*
90647c6ae99SBarry Smith          nc - number of components per grid point
90747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
90847c6ae99SBarry Smith 
90947c6ae99SBarry Smith   */
9101321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
91147c6ae99SBarry Smith   col = 2*s + 1;
912aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
913aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
91447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
91547c6ae99SBarry Smith 
91647c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
9171411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
9181411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
91947c6ae99SBarry Smith 
92047c6ae99SBarry Smith   /* determine the matrix preallocation information */
92147c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
92247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
92347c6ae99SBarry Smith 
9241321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9251321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
92647c6ae99SBarry Smith 
92747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
92847c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
92947c6ae99SBarry Smith 
9301321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9311321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
93247c6ae99SBarry Smith 
93347c6ae99SBarry Smith       for (k=0; k<nc; k++) {
93447c6ae99SBarry Smith         cnt  = 0;
93547c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
93647c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
93747c6ae99SBarry Smith             if (l || p) {
938aa219208SBarry Smith 	      if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
93947c6ae99SBarry Smith                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
94047c6ae99SBarry Smith 		  cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
94147c6ae99SBarry Smith 	      }
94247c6ae99SBarry Smith             } else {
94347c6ae99SBarry Smith 	      if (dfill) {
94447c6ae99SBarry Smith 		for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
94547c6ae99SBarry Smith 		  cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
94647c6ae99SBarry Smith 	      } else {
94747c6ae99SBarry Smith 		for (ifill_col=0; ifill_col<nc; ifill_col++)
94847c6ae99SBarry Smith 		  cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
94947c6ae99SBarry Smith 	      }
95047c6ae99SBarry Smith             }
95147c6ae99SBarry Smith 	  }
95247c6ae99SBarry Smith 	}
95347c6ae99SBarry Smith 	row = k + nc*(slot);
954784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
95547c6ae99SBarry Smith       }
95647c6ae99SBarry Smith     }
95747c6ae99SBarry Smith   }
95847c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
95947c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
96047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
961784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
962784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
96347c6ae99SBarry Smith 
96447c6ae99SBarry Smith   /*
96547c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
96647c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
96747c6ae99SBarry Smith     PETSc ordering.
96847c6ae99SBarry Smith   */
969fcfd50ebSBarry Smith   if (!da->prealloc_only) {
97047c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
97147c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
97247c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
97347c6ae99SBarry Smith 
9741321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9751321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
97647c6ae99SBarry Smith 
97747c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
97847c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
97947c6ae99SBarry Smith 
9801321219cSEthan Coon 	lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9811321219cSEthan Coon 	lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
98247c6ae99SBarry Smith 
98347c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
98447c6ae99SBarry Smith 	  cnt  = 0;
98547c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
98647c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
98747c6ae99SBarry Smith 	      if (l || p) {
988aa219208SBarry Smith 		if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
98947c6ae99SBarry Smith 		  for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
99047c6ae99SBarry Smith 		    cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
99147c6ae99SBarry Smith 		}
99247c6ae99SBarry Smith 	      } else {
99347c6ae99SBarry Smith 		if (dfill) {
99447c6ae99SBarry Smith 		  for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
99547c6ae99SBarry Smith 		    cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
99647c6ae99SBarry Smith 		} else {
99747c6ae99SBarry Smith 		  for (ifill_col=0; ifill_col<nc; ifill_col++)
99847c6ae99SBarry Smith 		    cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
99947c6ae99SBarry Smith 		}
100047c6ae99SBarry Smith 	      }
100147c6ae99SBarry Smith 	    }
100247c6ae99SBarry Smith 	  }
100347c6ae99SBarry Smith 	  row  = k + nc*(slot);
100447c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
100547c6ae99SBarry Smith 	}
100647c6ae99SBarry Smith       }
100747c6ae99SBarry Smith     }
100847c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
100947c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101047c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101147c6ae99SBarry Smith   }
101247c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
101347c6ae99SBarry Smith   PetscFunctionReturn(0);
101447c6ae99SBarry Smith }
101547c6ae99SBarry Smith 
101647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
101747c6ae99SBarry Smith 
101847c6ae99SBarry Smith #undef __FUNCT__
1019950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ"
1020950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM da,Mat J)
102147c6ae99SBarry Smith {
102247c6ae99SBarry Smith   PetscErrorCode         ierr;
102347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
102447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL;
102547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
102647c6ae99SBarry Smith   MPI_Comm               comm;
102747c6ae99SBarry Smith   PetscScalar            *values;
10281321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
102947c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1030aa219208SBarry Smith   DMDAStencilType        st;
103147c6ae99SBarry Smith 
103247c6ae99SBarry Smith   PetscFunctionBegin;
103347c6ae99SBarry Smith   /*
103447c6ae99SBarry Smith          nc - number of components per grid point
103547c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
103647c6ae99SBarry Smith 
103747c6ae99SBarry Smith   */
10381321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
103947c6ae99SBarry Smith   col    = 2*s + 1;
104047c6ae99SBarry Smith 
1041aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1042aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
104347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
104447c6ae99SBarry Smith 
104547c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
10461411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
10471411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
104847c6ae99SBarry Smith 
104947c6ae99SBarry Smith   /* determine the matrix preallocation information */
105047c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
105147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
10521321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10531321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
105447c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
10551321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10561321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
105747c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
10581321219cSEthan Coon 	kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
10591321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
106047c6ae99SBarry Smith 
106147c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
106247c6ae99SBarry Smith 
106347c6ae99SBarry Smith 	cnt  = 0;
106447c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
106547c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
106647c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
106747c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
1068aa219208SBarry Smith 		if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
106947c6ae99SBarry Smith 		  cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
107047c6ae99SBarry Smith 		}
107147c6ae99SBarry Smith 	      }
107247c6ae99SBarry Smith 	    }
107347c6ae99SBarry Smith 	  }
107447c6ae99SBarry Smith 	  rows[l] = l + nc*(slot);
107547c6ae99SBarry Smith 	}
1076784ac674SJed Brown 	ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
107747c6ae99SBarry Smith       }
107847c6ae99SBarry Smith     }
107947c6ae99SBarry Smith   }
1080f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
108147c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
108247c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
108347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1084784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1085784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
108647c6ae99SBarry Smith 
108747c6ae99SBarry Smith   /*
108847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
108947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
109047c6ae99SBarry Smith     PETSc ordering.
109147c6ae99SBarry Smith   */
1092fcfd50ebSBarry Smith   if (!da->prealloc_only) {
109347c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
109447c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
109547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
10961321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10971321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
109847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
10991321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
11001321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
110147c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
11021321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
11031321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
110447c6ae99SBarry Smith 
110547c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
110647c6ae99SBarry Smith 
110747c6ae99SBarry Smith 	  cnt  = 0;
110847c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
110947c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
111047c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
111147c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
1112aa219208SBarry Smith 		  if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
111347c6ae99SBarry Smith 		    cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
111447c6ae99SBarry Smith 		  }
111547c6ae99SBarry Smith 		}
111647c6ae99SBarry Smith 	      }
111747c6ae99SBarry Smith 	    }
111847c6ae99SBarry Smith 	    rows[l]      = l + nc*(slot);
111947c6ae99SBarry Smith 	  }
112047c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
112147c6ae99SBarry Smith 	}
112247c6ae99SBarry Smith       }
112347c6ae99SBarry Smith     }
112447c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
112547c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
112647c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
112747c6ae99SBarry Smith   }
112847c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
112947c6ae99SBarry Smith   PetscFunctionReturn(0);
113047c6ae99SBarry Smith }
113147c6ae99SBarry Smith 
113247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
113347c6ae99SBarry Smith 
113447c6ae99SBarry Smith #undef __FUNCT__
1135*ce308e1dSBarry Smith #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ_Fill"
1136*ce308e1dSBarry Smith PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM da,Mat J)
1137*ce308e1dSBarry Smith {
1138*ce308e1dSBarry Smith   PetscErrorCode         ierr;
1139*ce308e1dSBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1140*ce308e1dSBarry Smith   PetscInt               xs,nx,i,j,gxs,gnx,row,k,l;
1141*ce308e1dSBarry Smith   PetscInt               m,dim,s,*cols = PETSC_NULL,nc,col,cnt,*ocols;
1142*ce308e1dSBarry Smith   PetscInt               *ofill = dd->ofill;
1143*ce308e1dSBarry Smith   PetscScalar            *values;
1144*ce308e1dSBarry Smith   DMDABoundaryType       bx;
1145*ce308e1dSBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1146*ce308e1dSBarry Smith   PetscMPIInt            rank,size;
1147*ce308e1dSBarry Smith 
1148*ce308e1dSBarry Smith   PetscFunctionBegin;
1149*ce308e1dSBarry Smith   if (dd->bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"With fill provided not implemented with periodic boundary conditions");
1150*ce308e1dSBarry Smith   ierr = MPI_Comm_rank(((PetscObject)da)->comm,&rank);CHKERRQ(ierr);
1151*ce308e1dSBarry Smith   ierr = MPI_Comm_size(((PetscObject)da)->comm,&size);CHKERRQ(ierr);
1152*ce308e1dSBarry Smith 
1153*ce308e1dSBarry Smith   /*
1154*ce308e1dSBarry Smith          nc - number of components per grid point
1155*ce308e1dSBarry Smith          col - number of colors needed in one direction for single component problem
1156*ce308e1dSBarry Smith 
1157*ce308e1dSBarry Smith   */
1158*ce308e1dSBarry Smith   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
1159*ce308e1dSBarry Smith   col    = 2*s + 1;
1160*ce308e1dSBarry Smith 
1161*ce308e1dSBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1162*ce308e1dSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
1163*ce308e1dSBarry Smith 
1164*ce308e1dSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
1165*ce308e1dSBarry Smith   ierr = PetscMalloc2(nx*nc,PetscInt,&cols,nx*nc,PetscInt,&ocols);CHKERRQ(ierr);
1166*ce308e1dSBarry Smith   ierr = PetscMemzero(cols,nx*nc*sizeof(PetscInt));CHKERRQ(ierr);
1167*ce308e1dSBarry Smith   ierr = PetscMemzero(ocols,nx*nc*sizeof(PetscInt));CHKERRQ(ierr);
1168*ce308e1dSBarry Smith 
1169*ce308e1dSBarry Smith   /*
1170*ce308e1dSBarry Smith         note should be smaller for first and last process with no periodic
1171*ce308e1dSBarry Smith         does not handle dfill
1172*ce308e1dSBarry Smith   */
1173*ce308e1dSBarry Smith   cnt  = 0;
1174*ce308e1dSBarry Smith   /* coupling with process to the left */
1175*ce308e1dSBarry Smith   for (i=0; i<s; i++) {
1176*ce308e1dSBarry Smith     for (j=0; j<nc; j++){
1177*ce308e1dSBarry Smith       ocols[cnt] = ((!rank) ? 0 : (s - i)*(ofill[j+1] - ofill[j]));
1178*ce308e1dSBarry Smith       cols[cnt]  = nc + (s + i)*(ofill[j+1] - ofill[j]);
1179*ce308e1dSBarry Smith       cnt++;
1180*ce308e1dSBarry Smith     }
1181*ce308e1dSBarry Smith   }
1182*ce308e1dSBarry Smith   for (i=s; i<nx-s; i++) {
1183*ce308e1dSBarry Smith     for (j=0; j<nc; j++){
1184*ce308e1dSBarry Smith       cols[cnt]  = nc + 2*s*(ofill[j+1] - ofill[j]);
1185*ce308e1dSBarry Smith       cnt++;
1186*ce308e1dSBarry Smith     }
1187*ce308e1dSBarry Smith   }
1188*ce308e1dSBarry Smith  /* coupling with process to the right */
1189*ce308e1dSBarry Smith   for (i=nx-s; i<nx; i++){
1190*ce308e1dSBarry Smith     for (j=0; j<nc; j++){
1191*ce308e1dSBarry Smith       ocols[cnt] = ((rank == (size-1)) ? 0 : (i - nx + s + 1)*(ofill[j+1] - ofill[j]));
1192*ce308e1dSBarry Smith       cols[cnt]  = nc + (s + nx - i - 1)*(ofill[j+1] - ofill[j]);
1193*ce308e1dSBarry Smith       cnt++;
1194*ce308e1dSBarry Smith     }
1195*ce308e1dSBarry Smith   }
1196*ce308e1dSBarry Smith 
1197*ce308e1dSBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,cols);CHKERRQ(ierr);
1198*ce308e1dSBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,cols,0,ocols);CHKERRQ(ierr);
1199*ce308e1dSBarry Smith   ierr = PetscFree2(cols,ocols);CHKERRQ(ierr);
1200*ce308e1dSBarry Smith 
1201*ce308e1dSBarry Smith   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1202*ce308e1dSBarry Smith   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1203*ce308e1dSBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1204*ce308e1dSBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
1205*ce308e1dSBarry Smith 
1206*ce308e1dSBarry Smith   /*
1207*ce308e1dSBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1208*ce308e1dSBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1209*ce308e1dSBarry Smith     PETSc ordering.
1210*ce308e1dSBarry Smith   */
1211*ce308e1dSBarry Smith   if (!da->prealloc_only) {
1212*ce308e1dSBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1213*ce308e1dSBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1214*ce308e1dSBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1215*ce308e1dSBarry Smith 
1216*ce308e1dSBarry Smith     row  = xs*nc;
1217*ce308e1dSBarry Smith     /* coupling with process to the left */
1218*ce308e1dSBarry Smith     for (i=xs; i<xs+s; i++) {
1219*ce308e1dSBarry Smith       for (j=0; j<nc; j++){
1220*ce308e1dSBarry Smith         cnt = 0;
1221*ce308e1dSBarry Smith         if (rank) {
1222*ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1223*ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1224*ce308e1dSBarry Smith           }
1225*ce308e1dSBarry Smith         }
1226*ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1227*ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1228*ce308e1dSBarry Smith         }
1229*ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1230*ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1231*ce308e1dSBarry Smith         }
1232*ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1233*ce308e1dSBarry Smith         row++;
1234*ce308e1dSBarry Smith       }
1235*ce308e1dSBarry Smith     }
1236*ce308e1dSBarry Smith     for (i=xs+s; i<xs+nx-s; i++) {
1237*ce308e1dSBarry Smith       for (j=0; j<nc; j++){
1238*ce308e1dSBarry Smith         cnt = 0;
1239*ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1240*ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1241*ce308e1dSBarry Smith         }
1242*ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1243*ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1244*ce308e1dSBarry Smith         }
1245*ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1246*ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1247*ce308e1dSBarry Smith         }
1248*ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1249*ce308e1dSBarry Smith         row++;
1250*ce308e1dSBarry Smith       }
1251*ce308e1dSBarry Smith     }
1252*ce308e1dSBarry Smith     /* coupling with process to the right */
1253*ce308e1dSBarry Smith     for (i=xs+nx-s; i<xs+nx; i++){
1254*ce308e1dSBarry Smith       for (j=0; j<nc; j++){
1255*ce308e1dSBarry Smith         cnt = 0;
1256*ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1257*ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1258*ce308e1dSBarry Smith         }
1259*ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1260*ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1261*ce308e1dSBarry Smith         }
1262*ce308e1dSBarry Smith         if (rank < size-1) {
1263*ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1264*ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1265*ce308e1dSBarry Smith           }
1266*ce308e1dSBarry Smith         }
1267*ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1268*ce308e1dSBarry Smith         row++;
1269*ce308e1dSBarry Smith       }
1270*ce308e1dSBarry Smith     }
1271*ce308e1dSBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1272*ce308e1dSBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1273*ce308e1dSBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1274*ce308e1dSBarry Smith     ierr = PetscFree(cols);CHKERRQ(ierr);
1275*ce308e1dSBarry Smith   }
1276*ce308e1dSBarry Smith   PetscFunctionReturn(0);
1277*ce308e1dSBarry Smith }
1278*ce308e1dSBarry Smith 
1279*ce308e1dSBarry Smith /* ---------------------------------------------------------------------------------*/
1280*ce308e1dSBarry Smith 
1281*ce308e1dSBarry Smith #undef __FUNCT__
1282950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ"
1283950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM da,Mat J)
128447c6ae99SBarry Smith {
128547c6ae99SBarry Smith   PetscErrorCode         ierr;
128647c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
128747c6ae99SBarry Smith   PetscInt               m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l;
128847c6ae99SBarry Smith   PetscInt               istart,iend;
128947c6ae99SBarry Smith   PetscScalar            *values;
12901321219cSEthan Coon   DMDABoundaryType       bx;
129147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
129247c6ae99SBarry Smith 
129347c6ae99SBarry Smith   PetscFunctionBegin;
129447c6ae99SBarry Smith   /*
129547c6ae99SBarry Smith          nc - number of components per grid point
129647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
129747c6ae99SBarry Smith 
129847c6ae99SBarry Smith   */
12991321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
130047c6ae99SBarry Smith   col    = 2*s + 1;
130147c6ae99SBarry Smith 
1302aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1303aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
130447c6ae99SBarry Smith 
1305f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
130647c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
130747c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
130847c6ae99SBarry Smith 
13091411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
13101411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1311784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1312784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
131347c6ae99SBarry Smith 
131447c6ae99SBarry Smith   /*
131547c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
131647c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
131747c6ae99SBarry Smith     PETSc ordering.
131847c6ae99SBarry Smith   */
1319fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1320*ce308e1dSBarry Smith     ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
132147c6ae99SBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
132247c6ae99SBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
132347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
132447c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
132547c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
132647c6ae99SBarry Smith       slot   = i - gxs;
132747c6ae99SBarry Smith 
132847c6ae99SBarry Smith       cnt  = 0;
132947c6ae99SBarry Smith       for (l=0; l<nc; l++) {
133047c6ae99SBarry Smith 	for (i1=istart; i1<iend+1; i1++) {
133147c6ae99SBarry Smith 	  cols[cnt++] = l + nc*(slot + i1);
133247c6ae99SBarry Smith 	}
133347c6ae99SBarry Smith 	rows[l]      = l + nc*(slot);
133447c6ae99SBarry Smith       }
133547c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
133647c6ae99SBarry Smith     }
133747c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
133847c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
133947c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
134047c6ae99SBarry Smith     ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
1341*ce308e1dSBarry Smith   }
134247c6ae99SBarry Smith   PetscFunctionReturn(0);
134347c6ae99SBarry Smith }
134447c6ae99SBarry Smith 
134547c6ae99SBarry Smith #undef __FUNCT__
1346950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIBAIJ"
1347950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM da,Mat J)
134847c6ae99SBarry Smith {
134947c6ae99SBarry Smith   PetscErrorCode         ierr;
135047c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
135147c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
135247c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
135347c6ae99SBarry Smith   MPI_Comm               comm;
135447c6ae99SBarry Smith   PetscScalar            *values;
13551321219cSEthan Coon   DMDABoundaryType       bx,by;
1356aa219208SBarry Smith   DMDAStencilType        st;
135747c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
135847c6ae99SBarry Smith 
135947c6ae99SBarry Smith   PetscFunctionBegin;
136047c6ae99SBarry Smith   /*
136147c6ae99SBarry Smith      nc - number of components per grid point
136247c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
136347c6ae99SBarry Smith   */
13641321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
136547c6ae99SBarry Smith   col = 2*s + 1;
136647c6ae99SBarry Smith 
1367aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1368aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
136947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
137047c6ae99SBarry Smith 
137147c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
137247c6ae99SBarry Smith 
13731411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
13741411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
137547c6ae99SBarry Smith 
137647c6ae99SBarry Smith   /* determine the matrix preallocation information */
137747c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
137847c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
13791321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
13801321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
138147c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
13821321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
13831321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
138447c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
138547c6ae99SBarry Smith 
138647c6ae99SBarry Smith       /* Find block columns in block row */
138747c6ae99SBarry Smith       cnt  = 0;
138847c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
138947c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1390aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
139147c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
139247c6ae99SBarry Smith           }
139347c6ae99SBarry Smith         }
139447c6ae99SBarry Smith       }
1395784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
139647c6ae99SBarry Smith     }
139747c6ae99SBarry Smith   }
139847c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
139947c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
140047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
140147c6ae99SBarry Smith 
1402784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1403784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
140447c6ae99SBarry Smith 
140547c6ae99SBarry Smith   /*
140647c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
140747c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
140847c6ae99SBarry Smith     PETSc ordering.
140947c6ae99SBarry Smith   */
1410fcfd50ebSBarry Smith   if (!da->prealloc_only) {
141147c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
141247c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
141347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
14141321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14151321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
141647c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
14171321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14181321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
141947c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
142047c6ae99SBarry Smith 	cnt  = 0;
142147c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
142247c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1423aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
142447c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
142547c6ae99SBarry Smith             }
142647c6ae99SBarry Smith           }
142747c6ae99SBarry Smith         }
142847c6ae99SBarry Smith 	ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
142947c6ae99SBarry Smith       }
143047c6ae99SBarry Smith     }
143147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
143247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
143347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
143447c6ae99SBarry Smith   }
143547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
143647c6ae99SBarry Smith   PetscFunctionReturn(0);
143747c6ae99SBarry Smith }
143847c6ae99SBarry Smith 
143947c6ae99SBarry Smith #undef __FUNCT__
1440950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIBAIJ"
1441950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM da,Mat J)
144247c6ae99SBarry Smith {
144347c6ae99SBarry Smith   PetscErrorCode         ierr;
144447c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
144547c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
144647c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
144747c6ae99SBarry Smith   MPI_Comm               comm;
144847c6ae99SBarry Smith   PetscScalar            *values;
14491321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1450aa219208SBarry Smith   DMDAStencilType        st;
145147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
145247c6ae99SBarry Smith 
145347c6ae99SBarry Smith   PetscFunctionBegin;
145447c6ae99SBarry Smith   /*
145547c6ae99SBarry Smith          nc - number of components per grid point
145647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
145747c6ae99SBarry Smith 
145847c6ae99SBarry Smith   */
14591321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
146047c6ae99SBarry Smith   col    = 2*s + 1;
146147c6ae99SBarry Smith 
1462aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1463aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
146447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
146547c6ae99SBarry Smith 
146647c6ae99SBarry Smith   ierr  = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
146747c6ae99SBarry Smith 
14681411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
14691411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
147047c6ae99SBarry Smith 
147147c6ae99SBarry Smith   /* determine the matrix preallocation information */
147247c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
147347c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
14741321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14751321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
147647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
14771321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14781321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
147947c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
14801321219cSEthan Coon 	kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
14811321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
148247c6ae99SBarry Smith 
148347c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
148447c6ae99SBarry Smith 
148547c6ae99SBarry Smith 	/* Find block columns in block row */
148647c6ae99SBarry Smith 	cnt  = 0;
148747c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
148847c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
148947c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1490aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
149147c6ae99SBarry Smith 		cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
149247c6ae99SBarry Smith 	      }
149347c6ae99SBarry Smith 	    }
149447c6ae99SBarry Smith 	  }
149547c6ae99SBarry Smith 	}
1496784ac674SJed Brown 	ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
149747c6ae99SBarry Smith       }
149847c6ae99SBarry Smith     }
149947c6ae99SBarry Smith   }
150047c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
150147c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
150247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
150347c6ae99SBarry Smith 
1504784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1505784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
150647c6ae99SBarry Smith 
150747c6ae99SBarry Smith   /*
150847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
150947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
151047c6ae99SBarry Smith     PETSc ordering.
151147c6ae99SBarry Smith   */
1512fcfd50ebSBarry Smith   if (!da->prealloc_only) {
151347c6ae99SBarry Smith     ierr  = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
151447c6ae99SBarry Smith     ierr  = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
151547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
15161321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
15171321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
151847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
15191321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15201321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
152147c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
15221321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15231321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
152447c6ae99SBarry Smith 
152547c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
152647c6ae99SBarry Smith 
152747c6ae99SBarry Smith 	  cnt  = 0;
152847c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
152947c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
153047c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1531aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
153247c6ae99SBarry Smith                   cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
153347c6ae99SBarry Smith                 }
153447c6ae99SBarry Smith               }
153547c6ae99SBarry Smith             }
153647c6ae99SBarry Smith           }
153747c6ae99SBarry Smith 	  ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
153847c6ae99SBarry Smith 	}
153947c6ae99SBarry Smith       }
154047c6ae99SBarry Smith     }
154147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
154247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154447c6ae99SBarry Smith   }
154547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
154647c6ae99SBarry Smith   PetscFunctionReturn(0);
154747c6ae99SBarry Smith }
154847c6ae99SBarry Smith 
154947c6ae99SBarry Smith #undef __FUNCT__
155047c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular"
155147c6ae99SBarry Smith /*
155247c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
155347c6ae99SBarry Smith   identify in the local ordering with periodic domain.
155447c6ae99SBarry Smith */
155547c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
155647c6ae99SBarry Smith {
155747c6ae99SBarry Smith   PetscErrorCode ierr;
155847c6ae99SBarry Smith   PetscInt       i,n;
155947c6ae99SBarry Smith 
156047c6ae99SBarry Smith   PetscFunctionBegin;
156147c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr);
156247c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr);
156347c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
156447c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
156547c6ae99SBarry Smith   }
156647c6ae99SBarry Smith   *cnt = n;
156747c6ae99SBarry Smith   PetscFunctionReturn(0);
156847c6ae99SBarry Smith }
156947c6ae99SBarry Smith 
157047c6ae99SBarry Smith #undef __FUNCT__
1571950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPISBAIJ"
1572950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM da,Mat J)
157347c6ae99SBarry Smith {
157447c6ae99SBarry Smith   PetscErrorCode         ierr;
157547c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
157647c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
157747c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
157847c6ae99SBarry Smith   MPI_Comm               comm;
157947c6ae99SBarry Smith   PetscScalar            *values;
15801321219cSEthan Coon   DMDABoundaryType       bx,by;
1581aa219208SBarry Smith   DMDAStencilType        st;
158247c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
158347c6ae99SBarry Smith 
158447c6ae99SBarry Smith   PetscFunctionBegin;
158547c6ae99SBarry Smith   /*
158647c6ae99SBarry Smith      nc - number of components per grid point
158747c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
158847c6ae99SBarry Smith   */
15891321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
159047c6ae99SBarry Smith   col = 2*s + 1;
159147c6ae99SBarry Smith 
1592aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1593aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
159447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
159547c6ae99SBarry Smith 
159647c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
159747c6ae99SBarry Smith 
15981411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
15991411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
160047c6ae99SBarry Smith 
160147c6ae99SBarry Smith   /* determine the matrix preallocation information */
1602eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
160347c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
16041321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16051321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
160647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
16071321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16081321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
160947c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
161047c6ae99SBarry Smith 
161147c6ae99SBarry Smith       /* Find block columns in block row */
161247c6ae99SBarry Smith       cnt  = 0;
161347c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
161447c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1615aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) {
161647c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
161747c6ae99SBarry Smith           }
161847c6ae99SBarry Smith         }
161947c6ae99SBarry Smith       }
162047c6ae99SBarry Smith       ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
162147c6ae99SBarry Smith       ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
162247c6ae99SBarry Smith     }
162347c6ae99SBarry Smith   }
162447c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
162547c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
162647c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
162747c6ae99SBarry Smith 
1628784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1629784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
163047c6ae99SBarry Smith 
163147c6ae99SBarry Smith   /*
163247c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
163347c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
163447c6ae99SBarry Smith     PETSc ordering.
163547c6ae99SBarry Smith   */
1636fcfd50ebSBarry Smith   if (!da->prealloc_only) {
163747c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
163847c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
163947c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
16401321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16411321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
164247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
16431321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16441321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
164547c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
164647c6ae99SBarry Smith 
164747c6ae99SBarry Smith         /* Find block columns in block row */
164847c6ae99SBarry Smith         cnt  = 0;
164947c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
165047c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1651aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) {
165247c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
165347c6ae99SBarry Smith             }
165447c6ae99SBarry Smith           }
165547c6ae99SBarry Smith         }
165647c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
165747c6ae99SBarry Smith 	ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
165847c6ae99SBarry Smith       }
165947c6ae99SBarry Smith     }
166047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
166147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166347c6ae99SBarry Smith   }
166447c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
166547c6ae99SBarry Smith   PetscFunctionReturn(0);
166647c6ae99SBarry Smith }
166747c6ae99SBarry Smith 
166847c6ae99SBarry Smith #undef __FUNCT__
1669950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPISBAIJ"
1670950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM da,Mat J)
167147c6ae99SBarry Smith {
167247c6ae99SBarry Smith   PetscErrorCode         ierr;
167347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
167447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
167547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
167647c6ae99SBarry Smith   MPI_Comm               comm;
167747c6ae99SBarry Smith   PetscScalar            *values;
16781321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1679aa219208SBarry Smith   DMDAStencilType        st;
168047c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
168147c6ae99SBarry Smith 
168247c6ae99SBarry Smith   PetscFunctionBegin;
168347c6ae99SBarry Smith   /*
168447c6ae99SBarry Smith      nc - number of components per grid point
168547c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
168647c6ae99SBarry Smith   */
16871321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
168847c6ae99SBarry Smith   col = 2*s + 1;
168947c6ae99SBarry Smith 
1690aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1691aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
169247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
169347c6ae99SBarry Smith 
169447c6ae99SBarry Smith   /* create the matrix */
169547c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
169647c6ae99SBarry Smith 
16971411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
16981411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
169947c6ae99SBarry Smith 
170047c6ae99SBarry Smith   /* determine the matrix preallocation information */
1701eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
170247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
17031321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
17041321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
170547c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
17061321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
17071321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
170847c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
17091321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
17101321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
171147c6ae99SBarry Smith 
171247c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
171347c6ae99SBarry Smith 
171447c6ae99SBarry Smith 	/* Find block columns in block row */
171547c6ae99SBarry Smith 	cnt  = 0;
171647c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
171747c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
171847c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1719aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
172047c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
172147c6ae99SBarry Smith               }
172247c6ae99SBarry Smith             }
172347c6ae99SBarry Smith           }
172447c6ae99SBarry Smith         }
172547c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
172647c6ae99SBarry Smith         ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
172747c6ae99SBarry Smith       }
172847c6ae99SBarry Smith     }
172947c6ae99SBarry Smith   }
173047c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
173147c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
173247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
173347c6ae99SBarry Smith 
1734784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1735784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
173647c6ae99SBarry Smith 
173747c6ae99SBarry Smith   /*
173847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
173947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
174047c6ae99SBarry Smith     PETSc ordering.
174147c6ae99SBarry Smith   */
1742fcfd50ebSBarry Smith   if (!da->prealloc_only) {
174347c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
174447c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
174547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
17461321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
17471321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
174847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
17491321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
17501321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
175147c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
17521321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
17531321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
175447c6ae99SBarry Smith 
175547c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
175647c6ae99SBarry Smith 
175747c6ae99SBarry Smith 	  cnt  = 0;
175847c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
175947c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
176047c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1761aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
176247c6ae99SBarry Smith 		  cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
176347c6ae99SBarry Smith 		}
176447c6ae99SBarry Smith 	      }
176547c6ae99SBarry Smith 	    }
176647c6ae99SBarry Smith 	  }
176747c6ae99SBarry Smith           ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
176847c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
176947c6ae99SBarry Smith 	}
177047c6ae99SBarry Smith       }
177147c6ae99SBarry Smith     }
177247c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
177347c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177447c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177547c6ae99SBarry Smith   }
177647c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
177747c6ae99SBarry Smith   PetscFunctionReturn(0);
177847c6ae99SBarry Smith }
177947c6ae99SBarry Smith 
178047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
178147c6ae99SBarry Smith 
178247c6ae99SBarry Smith #undef __FUNCT__
1783950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ_Fill"
1784950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J)
178547c6ae99SBarry Smith {
178647c6ae99SBarry Smith   PetscErrorCode         ierr;
178747c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
178847c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz;
178947c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
179047c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
179147c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
179247c6ae99SBarry Smith   MPI_Comm               comm;
179347c6ae99SBarry Smith   PetscScalar            *values;
17941321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
179547c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1796aa219208SBarry Smith   DMDAStencilType        st;
179747c6ae99SBarry Smith 
179847c6ae99SBarry Smith   PetscFunctionBegin;
179947c6ae99SBarry Smith   /*
180047c6ae99SBarry Smith          nc - number of components per grid point
180147c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
180247c6ae99SBarry Smith 
180347c6ae99SBarry Smith   */
18041321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
180547c6ae99SBarry Smith   col    = 2*s + 1;
180666a15934SBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
180747c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
180866a15934SBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
180947c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
181066a15934SBarry Smith   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
181147c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
181247c6ae99SBarry Smith 
1813aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1814aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
181547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
181647c6ae99SBarry Smith 
181747c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
18181411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
18191411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
182047c6ae99SBarry Smith 
182147c6ae99SBarry Smith   /* determine the matrix preallocation information */
182247c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
182347c6ae99SBarry Smith 
182447c6ae99SBarry Smith 
182547c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
18261321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
18271321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
182847c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
18291321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
18301321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
183147c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
18321321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
18331321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
183447c6ae99SBarry Smith 
183547c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
183647c6ae99SBarry Smith 
183747c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
183847c6ae99SBarry Smith 	  cnt  = 0;
183947c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
184047c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
184147c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
184247c6ae99SBarry Smith 		if (ii || jj || kk) {
1843aa219208SBarry Smith 		  if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
184447c6ae99SBarry Smith 		    for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
184547c6ae99SBarry Smith 		      cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
184647c6ae99SBarry Smith 		  }
184747c6ae99SBarry Smith 		} else {
184847c6ae99SBarry Smith 		  if (dfill) {
184947c6ae99SBarry Smith 		    for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
185047c6ae99SBarry Smith 		      cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
185147c6ae99SBarry Smith 		  } else {
185247c6ae99SBarry Smith 		    for (ifill_col=0; ifill_col<nc; ifill_col++)
185347c6ae99SBarry Smith 		      cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
185447c6ae99SBarry Smith 		  }
185547c6ae99SBarry Smith 		}
185647c6ae99SBarry Smith 	      }
185747c6ae99SBarry Smith 	    }
185847c6ae99SBarry Smith 	  }
185947c6ae99SBarry Smith 	  row  = l + nc*(slot);
1860784ac674SJed Brown 	  ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
186147c6ae99SBarry Smith 	}
186247c6ae99SBarry Smith       }
186347c6ae99SBarry Smith     }
186447c6ae99SBarry Smith   }
186547c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
186647c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
186747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1868784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1869784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
187047c6ae99SBarry Smith 
187147c6ae99SBarry Smith   /*
187247c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
187347c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
187447c6ae99SBarry Smith     PETSc ordering.
187547c6ae99SBarry Smith   */
1876fcfd50ebSBarry Smith   if (!da->prealloc_only) {
187747c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
187847c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
187947c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
18801321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
18811321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
188247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
18831321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
18841321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
188547c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
18861321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
18871321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
188847c6ae99SBarry Smith 
188947c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
189047c6ae99SBarry Smith 
189147c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
189247c6ae99SBarry Smith 	    cnt  = 0;
189347c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
189447c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
189547c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
189647c6ae99SBarry Smith 		  if (ii || jj || kk) {
1897aa219208SBarry Smith 		    if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
189847c6ae99SBarry Smith 		      for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
189947c6ae99SBarry Smith 			cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
190047c6ae99SBarry Smith 		    }
190147c6ae99SBarry Smith 		  } else {
190247c6ae99SBarry Smith 		    if (dfill) {
190347c6ae99SBarry Smith 		      for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
190447c6ae99SBarry Smith 			cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
190547c6ae99SBarry Smith 		    } else {
190647c6ae99SBarry Smith 		      for (ifill_col=0; ifill_col<nc; ifill_col++)
190747c6ae99SBarry Smith 			cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
190847c6ae99SBarry Smith 		    }
190947c6ae99SBarry Smith 		  }
191047c6ae99SBarry Smith 		}
191147c6ae99SBarry Smith 	      }
191247c6ae99SBarry Smith 	    }
191347c6ae99SBarry Smith 	    row  = l + nc*(slot);
191447c6ae99SBarry Smith 	    ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
191547c6ae99SBarry Smith 	  }
191647c6ae99SBarry Smith 	}
191747c6ae99SBarry Smith       }
191847c6ae99SBarry Smith     }
191947c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
192047c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
192147c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
192247c6ae99SBarry Smith   }
192347c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
192447c6ae99SBarry Smith   PetscFunctionReturn(0);
192547c6ae99SBarry Smith }
1926