xref: /petsc/src/dm/impls/da/fdda.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
147c6ae99SBarry Smith 
2b45d2f2cSJed Brown #include <petsc-private/daimpl.h> /*I      "petscdmda.h"     I*/
307475bc1SBarry Smith #include <petscmat.h>
4b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
547c6ae99SBarry Smith 
6e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring*);
7e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring*);
8e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring*);
9e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM,ISColoringType,ISColoring*);
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith /*
1247c6ae99SBarry Smith    For ghost i that may be negative or greater than the upper bound this
1347c6ae99SBarry Smith   maps it into the 0:m-1 range using periodicity
1447c6ae99SBarry Smith */
1547c6ae99SBarry Smith #define SetInRange(i,m) ((i < 0) ? m+i : ((i >= m) ? i-m : i))
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith #undef __FUNCT__
18aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills_Private"
19ce308e1dSBarry Smith static PetscErrorCode DMDASetBlockFills_Private(const PetscInt *dfill,PetscInt w,PetscInt **rfill)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   PetscErrorCode ierr;
2247c6ae99SBarry Smith   PetscInt       i,j,nz,*fill;
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   PetscFunctionBegin;
2547c6ae99SBarry Smith   if (!dfill) PetscFunctionReturn(0);
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   /* count number nonzeros */
2847c6ae99SBarry Smith   nz = 0;
2947c6ae99SBarry Smith   for (i=0; i<w; i++) {
3047c6ae99SBarry Smith     for (j=0; j<w; j++) {
3147c6ae99SBarry Smith       if (dfill[w*i+j]) nz++;
3247c6ae99SBarry Smith     }
3347c6ae99SBarry Smith   }
3447c6ae99SBarry Smith   ierr = PetscMalloc((nz + w + 1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
3547c6ae99SBarry Smith   /* construct modified CSR storage of nonzero structure */
36ce308e1dSBarry Smith   /*  fill[0 -- w] marks starts of each row of column indices (and end of last row)
37ce308e1dSBarry Smith    so fill[1] - fill[0] gives number of nonzeros in first row etc */
3847c6ae99SBarry Smith   nz = w + 1;
3947c6ae99SBarry Smith   for (i=0; i<w; i++) {
4047c6ae99SBarry Smith     fill[i] = nz;
4147c6ae99SBarry Smith     for (j=0; j<w; j++) {
4247c6ae99SBarry Smith       if (dfill[w*i+j]) {
4347c6ae99SBarry Smith         fill[nz] = j;
4447c6ae99SBarry Smith         nz++;
4547c6ae99SBarry Smith       }
4647c6ae99SBarry Smith     }
4747c6ae99SBarry Smith   }
4847c6ae99SBarry Smith   fill[w] = nz;
4947c6ae99SBarry Smith 
5047c6ae99SBarry Smith   *rfill = fill;
5147c6ae99SBarry Smith   PetscFunctionReturn(0);
5247c6ae99SBarry Smith }
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith #undef __FUNCT__
55aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills"
5647c6ae99SBarry Smith /*@
57aa219208SBarry Smith     DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem
58950540a4SJed Brown     of the matrix returned by DMCreateMatrix().
5947c6ae99SBarry Smith 
60aa219208SBarry Smith     Logically Collective on DMDA
6147c6ae99SBarry Smith 
6247c6ae99SBarry Smith     Input Parameter:
6347c6ae99SBarry Smith +   da - the distributed array
640298fd71SBarry Smith .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
6547c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
6647c6ae99SBarry Smith 
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith     Level: developer
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     Notes: This only makes sense when you are doing multicomponent problems but using the
7147c6ae99SBarry Smith        MPIAIJ matrix format
7247c6ae99SBarry Smith 
7347c6ae99SBarry Smith            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
7447c6ae99SBarry Smith        representing coupling and 0 entries for missing coupling. For example
7547c6ae99SBarry Smith $             dfill[9] = {1, 0, 0,
7647c6ae99SBarry Smith $                         1, 1, 0,
7747c6ae99SBarry Smith $                         0, 1, 1}
7847c6ae99SBarry Smith        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
7947c6ae99SBarry Smith        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
8047c6ae99SBarry Smith        diagonal block).
8147c6ae99SBarry Smith 
82aa219208SBarry Smith      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
8347c6ae99SBarry Smith      can be represented in the dfill, ofill format
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith    Contributed by Glenn Hammond
8647c6ae99SBarry Smith 
878ddb5d8bSBarry Smith .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
8847c6ae99SBarry Smith 
8947c6ae99SBarry Smith @*/
90ce308e1dSBarry Smith PetscErrorCode  DMDASetBlockFills(DM da,const PetscInt *dfill,const PetscInt *ofill)
9147c6ae99SBarry Smith {
9247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
9347c6ae99SBarry Smith   PetscErrorCode ierr;
94ae4f298aSBarry Smith   PetscInt       i,k,cnt = 1;
9547c6ae99SBarry Smith 
9647c6ae99SBarry Smith   PetscFunctionBegin;
97aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr);
98aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr);
99ae4f298aSBarry Smith 
100ae4f298aSBarry Smith   /* ofillcount tracks the columns of ofill that have any nonzero in thems; the value in each location is the number of
101ae4f298aSBarry Smith    columns to the left with any nonzeros in them plus 1 */
102ae4f298aSBarry Smith   ierr = PetscMalloc(dd->w*sizeof(PetscBool),&dd->ofillcols);CHKERRQ(ierr);
103ae4f298aSBarry Smith   ierr = PetscMemzero(dd->ofillcols,dd->w*sizeof(PetscBool));CHKERRQ(ierr);
104ae4f298aSBarry Smith   for (i=0; i<dd->w; i++) {
105ae4f298aSBarry Smith     for (k=dd->ofill[i]; k<dd->ofill[i+1]; k++) dd->ofillcols[dd->ofill[k]] = 1;
106ae4f298aSBarry Smith   }
107ae4f298aSBarry Smith   for (i=0; i<dd->w; i++) {
108ae4f298aSBarry Smith     if (dd->ofillcols[i]) {
109ae4f298aSBarry Smith       dd->ofillcols[i] = cnt++;
110ae4f298aSBarry Smith     }
111ae4f298aSBarry Smith   }
11247c6ae99SBarry Smith   PetscFunctionReturn(0);
11347c6ae99SBarry Smith }
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith 
11647c6ae99SBarry Smith #undef __FUNCT__
117e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA"
11819fd82e9SBarry Smith PetscErrorCode  DMCreateColoring_DA(DM da,ISColoringType ctype,MatType mtype,ISColoring *coloring)
11947c6ae99SBarry Smith {
12047c6ae99SBarry Smith   PetscErrorCode   ierr;
12147c6ae99SBarry Smith   PetscInt         dim,m,n,p,nc;
1221321219cSEthan Coon   DMDABoundaryType bx,by,bz;
12347c6ae99SBarry Smith   MPI_Comm         comm;
12447c6ae99SBarry Smith   PetscMPIInt      size;
12547c6ae99SBarry Smith   PetscBool        isBAIJ;
12647c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith   PetscFunctionBegin;
12947c6ae99SBarry Smith   /*
13047c6ae99SBarry Smith                                   m
13147c6ae99SBarry Smith           ------------------------------------------------------
13247c6ae99SBarry Smith          |                                                     |
13347c6ae99SBarry Smith          |                                                     |
13447c6ae99SBarry Smith          |               ----------------------                |
13547c6ae99SBarry Smith          |               |                    |                |
13647c6ae99SBarry Smith       n  |           yn  |                    |                |
13747c6ae99SBarry Smith          |               |                    |                |
13847c6ae99SBarry Smith          |               .---------------------                |
13947c6ae99SBarry Smith          |             (xs,ys)     xn                          |
14047c6ae99SBarry Smith          |            .                                        |
14147c6ae99SBarry Smith          |         (gxs,gys)                                   |
14247c6ae99SBarry Smith          |                                                     |
14347c6ae99SBarry Smith           -----------------------------------------------------
14447c6ae99SBarry Smith   */
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith   /*
14747c6ae99SBarry Smith          nc - number of components per grid point
14847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
14947c6ae99SBarry Smith 
15047c6ae99SBarry Smith   */
1511321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&bx,&by,&bz,0);CHKERRQ(ierr);
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
15447c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
15547c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED) {
15647c6ae99SBarry Smith     if (size == 1) {
15747c6ae99SBarry Smith       ctype = IS_COLORING_GLOBAL;
15847c6ae99SBarry Smith     } else if (dim > 1) {
1591321219cSEthan Coon       if ((m==1 && bx == DMDA_BOUNDARY_PERIODIC) || (n==1 && by == DMDA_BOUNDARY_PERIODIC) || (p==1 && bz == DMDA_BOUNDARY_PERIODIC)) {
160*ce94432eSBarry Smith         SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"IS_COLORING_GHOSTED cannot be used for periodic boundary condition having both ends of the domain  on the same process");
16147c6ae99SBarry Smith       }
16247c6ae99SBarry Smith     }
16347c6ae99SBarry Smith   }
16447c6ae99SBarry Smith 
165aa219208SBarry Smith   /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ
16647c6ae99SBarry Smith      matrices is for the blocks, not the individual matrix elements  */
1674833614aSSatish Balay   ierr = PetscStrcmp(mtype,MATBAIJ,&isBAIJ);CHKERRQ(ierr);
1684833614aSSatish Balay   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);}
16947c6ae99SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);}
17047c6ae99SBarry Smith   if (isBAIJ) {
17147c6ae99SBarry Smith     dd->w  = 1;
17247c6ae99SBarry Smith     dd->xs = dd->xs/nc;
17347c6ae99SBarry Smith     dd->xe = dd->xe/nc;
17447c6ae99SBarry Smith     dd->Xs = dd->Xs/nc;
17547c6ae99SBarry Smith     dd->Xe = dd->Xe/nc;
17647c6ae99SBarry Smith   }
17747c6ae99SBarry Smith 
17847c6ae99SBarry Smith   /*
179aa219208SBarry Smith      We do not provide a getcoloring function in the DMDA operations because
180aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
18147c6ae99SBarry Smith    more low-level then matrices.
18247c6ae99SBarry Smith   */
18347c6ae99SBarry Smith   if (dim == 1) {
184e727c939SJed Brown     ierr = DMCreateColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
18547c6ae99SBarry Smith   } else if (dim == 2) {
186e727c939SJed Brown     ierr =  DMCreateColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
18747c6ae99SBarry Smith   } else if (dim == 3) {
188e727c939SJed Brown     ierr =  DMCreateColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
189*ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim);
19047c6ae99SBarry Smith   if (isBAIJ) {
19147c6ae99SBarry Smith     dd->w  = nc;
19247c6ae99SBarry Smith     dd->xs = dd->xs*nc;
19347c6ae99SBarry Smith     dd->xe = dd->xe*nc;
19447c6ae99SBarry Smith     dd->Xs = dd->Xs*nc;
19547c6ae99SBarry Smith     dd->Xe = dd->Xe*nc;
19647c6ae99SBarry Smith   }
19747c6ae99SBarry Smith   PetscFunctionReturn(0);
19847c6ae99SBarry Smith }
19947c6ae99SBarry Smith 
20047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith #undef __FUNCT__
203e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_MPIAIJ"
204e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
20547c6ae99SBarry Smith {
20647c6ae99SBarry Smith   PetscErrorCode   ierr;
20747c6ae99SBarry Smith   PetscInt         xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col;
20847c6ae99SBarry Smith   PetscInt         ncolors;
20947c6ae99SBarry Smith   MPI_Comm         comm;
2101321219cSEthan Coon   DMDABoundaryType bx,by;
211aa219208SBarry Smith   DMDAStencilType  st;
21247c6ae99SBarry Smith   ISColoringValue  *colors;
21347c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith   PetscFunctionBegin;
21647c6ae99SBarry Smith   /*
21747c6ae99SBarry Smith          nc - number of components per grid point
21847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
21947c6ae99SBarry Smith 
22047c6ae99SBarry Smith   */
2211321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
22247c6ae99SBarry Smith   col  = 2*s + 1;
223aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
224aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
22547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
22647c6ae99SBarry Smith 
22747c6ae99SBarry Smith   /* special case as taught to us by Paul Hovland */
228aa219208SBarry Smith   if (st == DMDA_STENCIL_STAR && s == 1) {
229e727c939SJed Brown     ierr = DMCreateColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
23047c6ae99SBarry Smith   } else {
23147c6ae99SBarry Smith 
232*ce94432eSBarry Smith     if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\
23347c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", m, col);
234*ce94432eSBarry Smith     if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\
23547c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", n, col);
23647c6ae99SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
23747c6ae99SBarry Smith       if (!dd->localcoloring) {
23847c6ae99SBarry Smith         ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
23947c6ae99SBarry Smith         ii   = 0;
24047c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
24147c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
24247c6ae99SBarry Smith             for (k=0; k<nc; k++) {
24347c6ae99SBarry Smith               colors[ii++] = k + nc*((i % col) + col*(j % col));
24447c6ae99SBarry Smith             }
24547c6ae99SBarry Smith           }
24647c6ae99SBarry Smith         }
24747c6ae99SBarry Smith         ncolors = nc + nc*(col-1 + col*(col-1));
24847c6ae99SBarry Smith         ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
24947c6ae99SBarry Smith       }
25047c6ae99SBarry Smith       *coloring = dd->localcoloring;
25147c6ae99SBarry Smith     } else if (ctype == IS_COLORING_GHOSTED) {
25247c6ae99SBarry Smith       if (!dd->ghostedcoloring) {
25347c6ae99SBarry Smith         ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
25447c6ae99SBarry Smith         ii   = 0;
25547c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
25647c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
25747c6ae99SBarry Smith             for (k=0; k<nc; k++) {
25847c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
25947c6ae99SBarry Smith               colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
26047c6ae99SBarry Smith             }
26147c6ae99SBarry Smith           }
26247c6ae99SBarry Smith         }
26347c6ae99SBarry Smith         ncolors = nc + nc*(col - 1 + col*(col-1));
26447c6ae99SBarry Smith         ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
26547c6ae99SBarry Smith         /* PetscIntView(ncolors,(PetscInt*)colors,0); */
26647c6ae99SBarry Smith 
26747c6ae99SBarry Smith         ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
26847c6ae99SBarry Smith       }
26947c6ae99SBarry Smith       *coloring = dd->ghostedcoloring;
270*ce94432eSBarry Smith     } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
27147c6ae99SBarry Smith   }
27247c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
27347c6ae99SBarry Smith   PetscFunctionReturn(0);
27447c6ae99SBarry Smith }
27547c6ae99SBarry Smith 
27647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
27747c6ae99SBarry Smith 
27847c6ae99SBarry Smith #undef __FUNCT__
279e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_3d_MPIAIJ"
280e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
28147c6ae99SBarry Smith {
28247c6ae99SBarry Smith   PetscErrorCode   ierr;
28347c6ae99SBarry Smith   PetscInt         xs,ys,nx,ny,i,j,gxs,gys,gnx,gny,m,n,p,dim,s,k,nc,col,zs,gzs,ii,l,nz,gnz,M,N,P;
28447c6ae99SBarry Smith   PetscInt         ncolors;
28547c6ae99SBarry Smith   MPI_Comm         comm;
2861321219cSEthan Coon   DMDABoundaryType bx,by,bz;
287aa219208SBarry Smith   DMDAStencilType  st;
28847c6ae99SBarry Smith   ISColoringValue  *colors;
28947c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
29047c6ae99SBarry Smith 
29147c6ae99SBarry Smith   PetscFunctionBegin;
29247c6ae99SBarry Smith   /*
29347c6ae99SBarry Smith          nc - number of components per grid point
29447c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
29547c6ae99SBarry Smith 
29647c6ae99SBarry Smith   */
2971321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
29847c6ae99SBarry Smith   col  = 2*s + 1;
299*ce94432eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
30047c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
301*ce94432eSBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
30247c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
303*ce94432eSBarry Smith   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
30447c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
30547c6ae99SBarry Smith 
306aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
307aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
30847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
30947c6ae99SBarry Smith 
31047c6ae99SBarry Smith   /* create the coloring */
31147c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
31247c6ae99SBarry Smith     if (!dd->localcoloring) {
31347c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
31447c6ae99SBarry Smith       ii   = 0;
31547c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
31647c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
31747c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
31847c6ae99SBarry Smith             for (l=0; l<nc; l++) {
31947c6ae99SBarry Smith               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
32047c6ae99SBarry Smith             }
32147c6ae99SBarry Smith           }
32247c6ae99SBarry Smith         }
32347c6ae99SBarry Smith       }
32447c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
32547c6ae99SBarry Smith       ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr);
32647c6ae99SBarry Smith     }
32747c6ae99SBarry Smith     *coloring = dd->localcoloring;
32847c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
32947c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
33047c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
33147c6ae99SBarry Smith       ii   = 0;
33247c6ae99SBarry Smith       for (k=gzs; k<gzs+gnz; k++) {
33347c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
33447c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
33547c6ae99SBarry Smith             for (l=0; l<nc; l++) {
33647c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
33747c6ae99SBarry Smith               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
33847c6ae99SBarry Smith             }
33947c6ae99SBarry Smith           }
34047c6ae99SBarry Smith         }
34147c6ae99SBarry Smith       }
34247c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
34347c6ae99SBarry Smith       ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
34447c6ae99SBarry Smith       ierr    = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
34547c6ae99SBarry Smith     }
34647c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
347*ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
34847c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
34947c6ae99SBarry Smith   PetscFunctionReturn(0);
35047c6ae99SBarry Smith }
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
35347c6ae99SBarry Smith 
35447c6ae99SBarry Smith #undef __FUNCT__
355e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_1d_MPIAIJ"
356e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
35747c6ae99SBarry Smith {
35847c6ae99SBarry Smith   PetscErrorCode   ierr;
35947c6ae99SBarry Smith   PetscInt         xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
36047c6ae99SBarry Smith   PetscInt         ncolors;
36147c6ae99SBarry Smith   MPI_Comm         comm;
3621321219cSEthan Coon   DMDABoundaryType bx;
36347c6ae99SBarry Smith   ISColoringValue  *colors;
36447c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
36547c6ae99SBarry Smith 
36647c6ae99SBarry Smith   PetscFunctionBegin;
36747c6ae99SBarry Smith   /*
36847c6ae99SBarry Smith          nc - number of components per grid point
36947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
37047c6ae99SBarry Smith 
37147c6ae99SBarry Smith   */
3721321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
37347c6ae99SBarry Smith   col  = 2*s + 1;
37447c6ae99SBarry Smith 
375*ce94432eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points %d is divisible\n\
37631e6f798SBarry Smith                                                           by 2*stencil_width + 1 %d\n",(int)m,(int)col);
37747c6ae99SBarry Smith 
378aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
379aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
38047c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
38147c6ae99SBarry Smith 
38247c6ae99SBarry Smith   /* create the coloring */
38347c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
38447c6ae99SBarry Smith     if (!dd->localcoloring) {
38547c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
386ae4f298aSBarry Smith       if (dd->ofillcols) {
387ae4f298aSBarry Smith         PetscInt tc = 0;
388ae4f298aSBarry Smith         for (i=0; i<nc; i++) tc += (PetscInt) (dd->ofillcols[i] > 0);
389ae4f298aSBarry Smith         i1 = 0;
390ae4f298aSBarry Smith         for (i=xs; i<xs+nx; i++) {
391ae4f298aSBarry Smith           for (l=0; l<nc; l++) {
392ae4f298aSBarry Smith             if (dd->ofillcols[l] && (i % col)) {
393ae4f298aSBarry Smith               colors[i1++] =  nc - 1 + tc*((i % col) - 1) + dd->ofillcols[l];
394ae4f298aSBarry Smith             } else {
395ae4f298aSBarry Smith               colors[i1++] = l;
396ae4f298aSBarry Smith             }
397ae4f298aSBarry Smith           }
398ae4f298aSBarry Smith         }
399ae4f298aSBarry Smith         ncolors = nc + 2*s*tc;
400ae4f298aSBarry Smith       } else {
40147c6ae99SBarry Smith         i1 = 0;
40247c6ae99SBarry Smith         for (i=xs; i<xs+nx; i++) {
40347c6ae99SBarry Smith           for (l=0; l<nc; l++) {
40447c6ae99SBarry Smith             colors[i1++] = l + nc*(i % col);
40547c6ae99SBarry Smith           }
40647c6ae99SBarry Smith         }
40747c6ae99SBarry Smith         ncolors = nc + nc*(col-1);
408ae4f298aSBarry Smith       }
40947c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr);
41047c6ae99SBarry Smith     }
41147c6ae99SBarry Smith     *coloring = dd->localcoloring;
41247c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
41347c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
41447c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
41547c6ae99SBarry Smith       i1   = 0;
41647c6ae99SBarry Smith       for (i=gxs; i<gxs+gnx; i++) {
41747c6ae99SBarry Smith         for (l=0; l<nc; l++) {
41847c6ae99SBarry Smith           /* the complicated stuff is to handle periodic boundaries */
41947c6ae99SBarry Smith           colors[i1++] = l + nc*(SetInRange(i,m) % col);
42047c6ae99SBarry Smith         }
42147c6ae99SBarry Smith       }
42247c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
42347c6ae99SBarry Smith       ierr    = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
42447c6ae99SBarry Smith       ierr    = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
42547c6ae99SBarry Smith     }
42647c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
427*ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
42847c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
42947c6ae99SBarry Smith   PetscFunctionReturn(0);
43047c6ae99SBarry Smith }
43147c6ae99SBarry Smith 
43247c6ae99SBarry Smith #undef __FUNCT__
433e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_5pt_MPIAIJ"
434e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
43547c6ae99SBarry Smith {
43647c6ae99SBarry Smith   PetscErrorCode   ierr;
43747c6ae99SBarry Smith   PetscInt         xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
43847c6ae99SBarry Smith   PetscInt         ncolors;
43947c6ae99SBarry Smith   MPI_Comm         comm;
4401321219cSEthan Coon   DMDABoundaryType bx,by;
44147c6ae99SBarry Smith   ISColoringValue  *colors;
44247c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
44347c6ae99SBarry Smith 
44447c6ae99SBarry Smith   PetscFunctionBegin;
44547c6ae99SBarry Smith   /*
44647c6ae99SBarry Smith          nc - number of components per grid point
44747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
44847c6ae99SBarry Smith 
44947c6ae99SBarry Smith   */
4501321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);CHKERRQ(ierr);
451aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
452aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
45347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
45447c6ae99SBarry Smith 
455*ce94432eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % 5)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible by 5\n");
456*ce94432eSBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % 5)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible by 5\n");
45747c6ae99SBarry Smith 
45847c6ae99SBarry Smith   /* create the coloring */
45947c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
46047c6ae99SBarry Smith     if (!dd->localcoloring) {
46147c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
46247c6ae99SBarry Smith       ii   = 0;
46347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
46447c6ae99SBarry Smith         for (i=xs; i<xs+nx; i++) {
46547c6ae99SBarry Smith           for (k=0; k<nc; k++) {
46647c6ae99SBarry Smith             colors[ii++] = k + nc*((3*j+i) % 5);
46747c6ae99SBarry Smith           }
46847c6ae99SBarry Smith         }
46947c6ae99SBarry Smith       }
47047c6ae99SBarry Smith       ncolors = 5*nc;
47147c6ae99SBarry Smith       ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
47247c6ae99SBarry Smith     }
47347c6ae99SBarry Smith     *coloring = dd->localcoloring;
47447c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
47547c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
47647c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
47747c6ae99SBarry Smith       ii = 0;
47847c6ae99SBarry Smith       for (j=gys; j<gys+gny; j++) {
47947c6ae99SBarry Smith         for (i=gxs; i<gxs+gnx; i++) {
48047c6ae99SBarry Smith           for (k=0; k<nc; k++) {
48147c6ae99SBarry Smith             colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
48247c6ae99SBarry Smith           }
48347c6ae99SBarry Smith         }
48447c6ae99SBarry Smith       }
48547c6ae99SBarry Smith       ncolors = 5*nc;
48647c6ae99SBarry Smith       ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
48747c6ae99SBarry Smith       ierr    = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
48847c6ae99SBarry Smith     }
48947c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
490*ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
49147c6ae99SBarry Smith   PetscFunctionReturn(0);
49247c6ae99SBarry Smith }
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith /* =========================================================================== */
495950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM,Mat);
496ce308e1dSBarry Smith extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM,Mat);
497950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM,Mat);
498950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM,Mat);
499950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM,Mat);
500950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM,Mat);
501950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM,Mat);
502950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM,Mat);
503950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM,Mat);
504950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM,Mat);
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith #undef __FUNCT__
507c688c046SMatthew G Knepley #define __FUNCT__ "MatSetupDM"
5088bbdbebaSMatthew G Knepley /*@C
509c688c046SMatthew G Knepley    MatSetupDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith    Logically Collective on Mat
51247c6ae99SBarry Smith 
51347c6ae99SBarry Smith    Input Parameters:
51447c6ae99SBarry Smith +  mat - the matrix
51547c6ae99SBarry Smith -  da - the da
51647c6ae99SBarry Smith 
51747c6ae99SBarry Smith    Level: intermediate
51847c6ae99SBarry Smith 
51947c6ae99SBarry Smith @*/
520c688c046SMatthew G Knepley PetscErrorCode MatSetupDM(Mat mat,DM da)
52147c6ae99SBarry Smith {
52247c6ae99SBarry Smith   PetscErrorCode ierr;
52347c6ae99SBarry Smith 
52447c6ae99SBarry Smith   PetscFunctionBegin;
52547c6ae99SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
52647c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
527c688c046SMatthew G Knepley   ierr = PetscTryMethod(mat,"MatSetupDM_C",(Mat,DM),(mat,da));CHKERRQ(ierr);
52847c6ae99SBarry Smith   PetscFunctionReturn(0);
52947c6ae99SBarry Smith }
53047c6ae99SBarry Smith 
53147c6ae99SBarry Smith EXTERN_C_BEGIN
53247c6ae99SBarry Smith #undef __FUNCT__
53347c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA"
5347087cfbeSBarry Smith PetscErrorCode  MatView_MPI_DA(Mat A,PetscViewer viewer)
53547c6ae99SBarry Smith {
5369a42bb27SBarry Smith   DM                da;
53747c6ae99SBarry Smith   PetscErrorCode    ierr;
53847c6ae99SBarry Smith   const char        *prefix;
53947c6ae99SBarry Smith   Mat               Anatural;
54047c6ae99SBarry Smith   AO                ao;
54147c6ae99SBarry Smith   PetscInt          rstart,rend,*petsc,i;
54247c6ae99SBarry Smith   IS                is;
54347c6ae99SBarry Smith   MPI_Comm          comm;
54474388724SJed Brown   PetscViewerFormat format;
54547c6ae99SBarry Smith 
54647c6ae99SBarry Smith   PetscFunctionBegin;
54774388724SJed Brown   /* Check whether we are just printing info, in which case MatView() already viewed everything we wanted to view */
54874388724SJed Brown   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
54974388724SJed Brown   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
55074388724SJed Brown 
55147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
552c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
553*ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
55447c6ae99SBarry Smith 
555aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
55647c6ae99SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
55747c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr);
55847c6ae99SBarry Smith   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
55947c6ae99SBarry Smith   ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr);
56047c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith   /* call viewer on natural ordering */
56347c6ae99SBarry Smith   ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr);
564fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
56547c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr);
56647c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr);
56747c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr);
56847c6ae99SBarry Smith   ierr = MatView(Anatural,viewer);CHKERRQ(ierr);
569fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
57047c6ae99SBarry Smith   PetscFunctionReturn(0);
57147c6ae99SBarry Smith }
57247c6ae99SBarry Smith EXTERN_C_END
57347c6ae99SBarry Smith 
57447c6ae99SBarry Smith EXTERN_C_BEGIN
57547c6ae99SBarry Smith #undef __FUNCT__
57647c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA"
5777087cfbeSBarry Smith PetscErrorCode  MatLoad_MPI_DA(Mat A,PetscViewer viewer)
57847c6ae99SBarry Smith {
5799a42bb27SBarry Smith   DM             da;
58047c6ae99SBarry Smith   PetscErrorCode ierr;
58147c6ae99SBarry Smith   Mat            Anatural,Aapp;
58247c6ae99SBarry Smith   AO             ao;
58347c6ae99SBarry Smith   PetscInt       rstart,rend,*app,i;
58447c6ae99SBarry Smith   IS             is;
58547c6ae99SBarry Smith   MPI_Comm       comm;
58647c6ae99SBarry Smith 
58747c6ae99SBarry Smith   PetscFunctionBegin;
58847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
589c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
590*ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
59147c6ae99SBarry Smith 
59247c6ae99SBarry Smith   /* Load the matrix in natural ordering */
593*ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&Anatural);CHKERRQ(ierr);
59447c6ae99SBarry Smith   ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr);
59547c6ae99SBarry Smith   ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
59647c6ae99SBarry Smith   ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr);
59747c6ae99SBarry Smith 
59847c6ae99SBarry Smith   /* Map natural ordering to application ordering and create IS */
599aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
60047c6ae99SBarry Smith   ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr);
60147c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr);
60247c6ae99SBarry Smith   for (i=rstart; i<rend; i++) app[i-rstart] = i;
60347c6ae99SBarry Smith   ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr);
60447c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
60547c6ae99SBarry Smith 
60647c6ae99SBarry Smith   /* Do permutation and replace header */
60747c6ae99SBarry Smith   ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr);
60847c6ae99SBarry Smith   ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr);
609fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
610fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
61147c6ae99SBarry Smith   PetscFunctionReturn(0);
61247c6ae99SBarry Smith }
61347c6ae99SBarry Smith EXTERN_C_END
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith #undef __FUNCT__
616950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA"
61719fd82e9SBarry Smith PetscErrorCode DMCreateMatrix_DA(DM da, MatType mtype,Mat *J)
61847c6ae99SBarry Smith {
61947c6ae99SBarry Smith   PetscErrorCode ierr;
62047c6ae99SBarry Smith   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
62147c6ae99SBarry Smith   Mat            A;
62247c6ae99SBarry Smith   MPI_Comm       comm;
62319fd82e9SBarry Smith   MatType        Atype;
62437d0c07bSMatthew G Knepley   PetscSection   section, sectionGlobal;
6250298fd71SBarry Smith   void           (*aij)(void)=NULL,(*baij)(void)=NULL,(*sbaij)(void)=NULL;
62647c6ae99SBarry Smith   MatType        ttype[256];
62747c6ae99SBarry Smith   PetscBool      flg;
62847c6ae99SBarry Smith   PetscMPIInt    size;
62947c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
63047c6ae99SBarry Smith 
63147c6ae99SBarry Smith   PetscFunctionBegin;
632519f805aSKarl Rupp #if !defined(PETSC_USE_DYNAMIC_LIBRARIES)
6330298fd71SBarry Smith   ierr = MatInitializePackage(NULL);CHKERRQ(ierr);
63447c6ae99SBarry Smith #endif
6355da5aae0SJed Brown   if (!mtype) mtype = MATAIJ;
63647c6ae99SBarry Smith   ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr);
637*ce94432eSBarry Smith   ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)da),((PetscObject)da)->prefix,"DMDA options","Mat");CHKERRQ(ierr);
638dd85299cSBarry Smith   ierr = PetscOptionsList("-dm_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr);
63947c6ae99SBarry Smith   ierr = PetscOptionsEnd();
64047c6ae99SBarry Smith 
64137d0c07bSMatthew G Knepley   ierr = DMGetDefaultSection(da, &section);CHKERRQ(ierr);
64237d0c07bSMatthew G Knepley   if (section) {
64337d0c07bSMatthew G Knepley     PetscInt  bs = -1;
64437d0c07bSMatthew G Knepley     PetscInt  localSize;
64537d0c07bSMatthew G Knepley     PetscBool isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isSymmetric;
64637d0c07bSMatthew G Knepley 
64737d0c07bSMatthew G Knepley     ierr = DMGetDefaultGlobalSection(da, &sectionGlobal);CHKERRQ(ierr);
64837d0c07bSMatthew G Knepley     ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize);CHKERRQ(ierr);
64937d0c07bSMatthew G Knepley     ierr = MatCreate(((PetscObject) da)->comm, J);CHKERRQ(ierr);
65037d0c07bSMatthew G Knepley     ierr = MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
65137d0c07bSMatthew G Knepley     ierr = MatSetType(*J, mtype);CHKERRQ(ierr);
65237d0c07bSMatthew G Knepley     ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
65337d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSHELL, &isShell);CHKERRQ(ierr);
65437d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATBAIJ, &isBlock);CHKERRQ(ierr);
65537d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSEQBAIJ, &isSeqBlock);CHKERRQ(ierr);
65637d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATMPIBAIJ, &isMPIBlock);CHKERRQ(ierr);
65737d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr);
65837d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr);
65937d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr);
66037d0c07bSMatthew G Knepley     /* Check for symmetric storage */
66137d0c07bSMatthew G Knepley     isSymmetric = (PetscBool) (isSymBlock || isSymSeqBlock || isSymMPIBlock);
66237d0c07bSMatthew G Knepley     if (isSymmetric) {
66337d0c07bSMatthew G Knepley       ierr = MatSetOption(*J, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);
66437d0c07bSMatthew G Knepley     }
66537d0c07bSMatthew G Knepley     if (!isShell) {
66637d0c07bSMatthew G Knepley       PetscInt *dnz, *onz, *dnzu, *onzu, bsLocal;
66737d0c07bSMatthew G Knepley 
66837d0c07bSMatthew G Knepley       if (bs < 0) {
66937d0c07bSMatthew G Knepley         if (isBlock || isSeqBlock || isMPIBlock || isSymBlock || isSymSeqBlock || isSymMPIBlock) {
67037d0c07bSMatthew G Knepley           PetscInt pStart, pEnd, p, dof;
67137d0c07bSMatthew G Knepley 
67237d0c07bSMatthew G Knepley           ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
67337d0c07bSMatthew G Knepley           for (p = pStart; p < pEnd; ++p) {
67437d0c07bSMatthew G Knepley             ierr = PetscSectionGetDof(sectionGlobal, p, &dof);CHKERRQ(ierr);
67537d0c07bSMatthew G Knepley             if (dof) {
67637d0c07bSMatthew G Knepley               bs = dof;
67737d0c07bSMatthew G Knepley               break;
67837d0c07bSMatthew G Knepley             }
67937d0c07bSMatthew G Knepley           }
68037d0c07bSMatthew G Knepley         } else {
68137d0c07bSMatthew G Knepley           bs = 1;
68237d0c07bSMatthew G Knepley         }
68337d0c07bSMatthew G Knepley         /* Must have same blocksize on all procs (some might have no points) */
68437d0c07bSMatthew G Knepley         bsLocal = bs;
68537d0c07bSMatthew G Knepley         ierr    = MPI_Allreduce(&bsLocal, &bs, 1, MPIU_INT, MPI_MAX, ((PetscObject) da)->comm);CHKERRQ(ierr);
68637d0c07bSMatthew G Knepley       }
68737d0c07bSMatthew G Knepley       ierr = PetscMalloc4(localSize/bs, PetscInt, &dnz, localSize/bs, PetscInt, &onz, localSize/bs, PetscInt, &dnzu, localSize/bs, PetscInt, &onzu);CHKERRQ(ierr);
68837d0c07bSMatthew G Knepley       ierr = PetscMemzero(dnz,  localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
68937d0c07bSMatthew G Knepley       ierr = PetscMemzero(onz,  localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
69037d0c07bSMatthew G Knepley       ierr = PetscMemzero(dnzu, localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
69137d0c07bSMatthew G Knepley       ierr = PetscMemzero(onzu, localSize/bs * sizeof(PetscInt));CHKERRQ(ierr);
692552f7358SJed Brown       /* ierr = DMPlexPreallocateOperator(dm, bs, section, sectionGlobal, dnz, onz, dnzu, onzu, *J, fillMatrix);CHKERRQ(ierr); */
69337d0c07bSMatthew G Knepley       ierr = PetscFree4(dnz, onz, dnzu, onzu);CHKERRQ(ierr);
69437d0c07bSMatthew G Knepley     }
69537d0c07bSMatthew G Knepley   }
69647c6ae99SBarry Smith   /*
69747c6ae99SBarry Smith                                   m
69847c6ae99SBarry Smith           ------------------------------------------------------
69947c6ae99SBarry Smith          |                                                     |
70047c6ae99SBarry Smith          |                                                     |
70147c6ae99SBarry Smith          |               ----------------------                |
70247c6ae99SBarry Smith          |               |                    |                |
70347c6ae99SBarry Smith       n  |           ny  |                    |                |
70447c6ae99SBarry Smith          |               |                    |                |
70547c6ae99SBarry Smith          |               .---------------------                |
70647c6ae99SBarry Smith          |             (xs,ys)     nx                          |
70747c6ae99SBarry Smith          |            .                                        |
70847c6ae99SBarry Smith          |         (gxs,gys)                                   |
70947c6ae99SBarry Smith          |                                                     |
71047c6ae99SBarry Smith           -----------------------------------------------------
71147c6ae99SBarry Smith   */
71247c6ae99SBarry Smith 
71347c6ae99SBarry Smith   /*
71447c6ae99SBarry Smith          nc - number of components per grid point
71547c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
71647c6ae99SBarry Smith 
71747c6ae99SBarry Smith   */
718e30e807fSPeter Brune   M   = dd->M;
719e30e807fSPeter Brune   N   = dd->N;
720e30e807fSPeter Brune   P   = dd->P;
721e30e807fSPeter Brune   dim = dd->dim;
722e30e807fSPeter Brune   dof = dd->w;
723e30e807fSPeter Brune   /* ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); */
724aa219208SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
72547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
72647c6ae99SBarry Smith   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
72747c6ae99SBarry Smith   ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr);
72819fd82e9SBarry Smith   ierr = MatSetType(A,(MatType)ttype);CHKERRQ(ierr);
72995ee5b0eSBarry Smith   ierr = MatSetDM(A,da);CHKERRQ(ierr);
73047c6ae99SBarry Smith   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
73147c6ae99SBarry Smith   ierr = MatGetType(A,&Atype);CHKERRQ(ierr);
73247c6ae99SBarry Smith   /*
733aa219208SBarry Smith      We do not provide a getmatrix function in the DMDA operations because
734aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
73547c6ae99SBarry Smith    more low-level than matrices. This is kind of cheating but, cause sometimes
736aa219208SBarry Smith    we think of DMDA has higher level than matrices.
73747c6ae99SBarry Smith 
73847c6ae99SBarry Smith      We could switch based on Atype (or mtype), but we do not since the
73947c6ae99SBarry Smith    specialized setting routines depend only the particular preallocation
74047c6ae99SBarry Smith    details of the matrix, not the type itself.
74147c6ae99SBarry Smith   */
74247c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
74347c6ae99SBarry Smith   if (!aij) {
74447c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
74547c6ae99SBarry Smith   }
74647c6ae99SBarry Smith   if (!aij) {
74747c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
74847c6ae99SBarry Smith     if (!baij) {
74947c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
75047c6ae99SBarry Smith     }
75147c6ae99SBarry Smith     if (!baij) {
75247c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
75347c6ae99SBarry Smith       if (!sbaij) {
75447c6ae99SBarry Smith         ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
75547c6ae99SBarry Smith       }
75647c6ae99SBarry Smith     }
75747c6ae99SBarry Smith   }
75847c6ae99SBarry Smith   if (aij) {
75947c6ae99SBarry Smith     if (dim == 1) {
760ce308e1dSBarry Smith       if (dd->ofill) {
761ce308e1dSBarry Smith         ierr = DMCreateMatrix_DA_1d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
762ce308e1dSBarry Smith       } else {
763950540a4SJed Brown         ierr = DMCreateMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr);
764ce308e1dSBarry Smith       }
76547c6ae99SBarry Smith     } else if (dim == 2) {
76647c6ae99SBarry Smith       if (dd->ofill) {
767950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
76847c6ae99SBarry Smith       } else {
769950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr);
77047c6ae99SBarry Smith       }
77147c6ae99SBarry Smith     } else if (dim == 3) {
77247c6ae99SBarry Smith       if (dd->ofill) {
773950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
77447c6ae99SBarry Smith       } else {
775950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr);
77647c6ae99SBarry Smith       }
77747c6ae99SBarry Smith     }
77847c6ae99SBarry Smith   } else if (baij) {
77947c6ae99SBarry Smith     if (dim == 2) {
780950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr);
78147c6ae99SBarry Smith     } else if (dim == 3) {
782950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr);
783*ce94432eSBarry Smith     } else SETERRQ3(PetscObjectComm((PetscObject)da),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);
78447c6ae99SBarry Smith   } else if (sbaij) {
78547c6ae99SBarry Smith     if (dim == 2) {
786950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr);
78747c6ae99SBarry Smith     } else if (dim == 3) {
788950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr);
789*ce94432eSBarry Smith     } else SETERRQ3(PetscObjectComm((PetscObject)da),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);
790869776cdSLisandro Dalcin   } else {
791869776cdSLisandro Dalcin     ISLocalToGlobalMapping ltog,ltogb;
792869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
793869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
7942949035bSJed Brown     ierr = MatSetUp(A);CHKERRQ(ierr);
795869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMapping(A,ltog,ltog);CHKERRQ(ierr);
796869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMappingBlock(A,ltogb,ltogb);CHKERRQ(ierr);
79747c6ae99SBarry Smith   }
798aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
79947c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
800c688c046SMatthew G Knepley   ierr = MatSetDM(A,da);CHKERRQ(ierr);
80147c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
80247c6ae99SBarry Smith   if (size > 1) {
80347c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
80447c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MPI_DA);CHKERRQ(ierr);
80547c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void))MatLoad_MPI_DA);CHKERRQ(ierr);
80647c6ae99SBarry Smith   }
80747c6ae99SBarry Smith   *J = A;
80847c6ae99SBarry Smith   PetscFunctionReturn(0);
80947c6ae99SBarry Smith }
81047c6ae99SBarry Smith 
81147c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
81247c6ae99SBarry Smith #undef __FUNCT__
813950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ"
814950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM da,Mat J)
81547c6ae99SBarry Smith {
81647c6ae99SBarry Smith   PetscErrorCode         ierr;
8170298fd71SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny,m,n,dim,s,*cols = NULL,k,nc,*rows = NULL,col,cnt,l,p;
81847c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
81947c6ae99SBarry Smith   MPI_Comm               comm;
82047c6ae99SBarry Smith   PetscScalar            *values;
8211321219cSEthan Coon   DMDABoundaryType       bx,by;
82247c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
823aa219208SBarry Smith   DMDAStencilType        st;
82447c6ae99SBarry Smith 
82547c6ae99SBarry Smith   PetscFunctionBegin;
82647c6ae99SBarry Smith   /*
82747c6ae99SBarry Smith          nc - number of components per grid point
82847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith   */
8311321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
83247c6ae99SBarry Smith   col  = 2*s + 1;
833aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
834aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
83547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
83647c6ae99SBarry Smith 
83747c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
8381411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
8391411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
84047c6ae99SBarry Smith 
84147c6ae99SBarry Smith   /* determine the matrix preallocation information */
84247c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
84347c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
84447c6ae99SBarry Smith 
8451321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8461321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
84947c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
85047c6ae99SBarry Smith 
8511321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8521321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith       cnt = 0;
85547c6ae99SBarry Smith       for (k=0; k<nc; k++) {
85647c6ae99SBarry Smith         for (l=lstart; l<lend+1; l++) {
85747c6ae99SBarry Smith           for (p=pstart; p<pend+1; p++) {
858aa219208SBarry Smith             if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
85947c6ae99SBarry Smith               cols[cnt++] = k + nc*(slot + gnx*l + p);
86047c6ae99SBarry Smith             }
86147c6ae99SBarry Smith           }
86247c6ae99SBarry Smith         }
86347c6ae99SBarry Smith         rows[k] = k + nc*(slot);
86447c6ae99SBarry Smith       }
865784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
86647c6ae99SBarry Smith     }
86747c6ae99SBarry Smith   }
868f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
86947c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
87047c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
87147c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
87247c6ae99SBarry Smith 
873784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
874784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
87547c6ae99SBarry Smith 
87647c6ae99SBarry Smith   /*
87747c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
87847c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
87947c6ae99SBarry Smith     PETSc ordering.
88047c6ae99SBarry Smith   */
881fcfd50ebSBarry Smith   if (!da->prealloc_only) {
88247c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
88347c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
88447c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
88547c6ae99SBarry Smith 
8861321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8871321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
88847c6ae99SBarry Smith 
88947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
89047c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
89147c6ae99SBarry Smith 
8921321219cSEthan Coon         lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8931321219cSEthan Coon         lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
89447c6ae99SBarry Smith 
89547c6ae99SBarry Smith         cnt = 0;
89647c6ae99SBarry Smith         for (k=0; k<nc; k++) {
89747c6ae99SBarry Smith           for (l=lstart; l<lend+1; l++) {
89847c6ae99SBarry Smith             for (p=pstart; p<pend+1; p++) {
899aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
90047c6ae99SBarry Smith                 cols[cnt++] = k + nc*(slot + gnx*l + p);
90147c6ae99SBarry Smith               }
90247c6ae99SBarry Smith             }
90347c6ae99SBarry Smith           }
90447c6ae99SBarry Smith           rows[k] = k + nc*(slot);
90547c6ae99SBarry Smith         }
90647c6ae99SBarry Smith         ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
90747c6ae99SBarry Smith       }
90847c6ae99SBarry Smith     }
90947c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
91047c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91147c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91247c6ae99SBarry Smith   }
91347c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
91447c6ae99SBarry Smith   PetscFunctionReturn(0);
91547c6ae99SBarry Smith }
91647c6ae99SBarry Smith 
91747c6ae99SBarry Smith #undef __FUNCT__
918950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ_Fill"
919950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J)
92047c6ae99SBarry Smith {
92147c6ae99SBarry Smith   PetscErrorCode         ierr;
92247c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
92347c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p;
92447c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
92547c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
92647c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
92747c6ae99SBarry Smith   MPI_Comm               comm;
92847c6ae99SBarry Smith   PetscScalar            *values;
9291321219cSEthan Coon   DMDABoundaryType       bx,by;
93047c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
931aa219208SBarry Smith   DMDAStencilType        st;
93247c6ae99SBarry Smith 
93347c6ae99SBarry Smith   PetscFunctionBegin;
93447c6ae99SBarry Smith   /*
93547c6ae99SBarry Smith          nc - number of components per grid point
93647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
93747c6ae99SBarry Smith 
93847c6ae99SBarry Smith   */
9391321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
94047c6ae99SBarry Smith   col  = 2*s + 1;
941aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
942aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
94347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
94447c6ae99SBarry Smith 
94547c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
9461411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
9471411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
94847c6ae99SBarry Smith 
94947c6ae99SBarry Smith   /* determine the matrix preallocation information */
95047c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
95147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
95247c6ae99SBarry Smith 
9531321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9541321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
95547c6ae99SBarry Smith 
95647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
95747c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
95847c6ae99SBarry Smith 
9591321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9601321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
96147c6ae99SBarry Smith 
96247c6ae99SBarry Smith       for (k=0; k<nc; k++) {
96347c6ae99SBarry Smith         cnt = 0;
96447c6ae99SBarry Smith         for (l=lstart; l<lend+1; l++) {
96547c6ae99SBarry Smith           for (p=pstart; p<pend+1; p++) {
96647c6ae99SBarry Smith             if (l || p) {
967aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
9688865f1eaSKarl Rupp                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p);
96947c6ae99SBarry Smith               }
97047c6ae99SBarry Smith             } else {
97147c6ae99SBarry Smith               if (dfill) {
9728865f1eaSKarl Rupp                 for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p);
97347c6ae99SBarry Smith               } else {
9748865f1eaSKarl Rupp                 for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + gnx*l + p);
97547c6ae99SBarry Smith               }
97647c6ae99SBarry Smith             }
97747c6ae99SBarry Smith           }
97847c6ae99SBarry Smith         }
97947c6ae99SBarry Smith         row  = k + nc*(slot);
980784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
98147c6ae99SBarry Smith       }
98247c6ae99SBarry Smith     }
98347c6ae99SBarry Smith   }
98447c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
98547c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
98647c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
987784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
988784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
98947c6ae99SBarry Smith 
99047c6ae99SBarry Smith   /*
99147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
99247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
99347c6ae99SBarry Smith     PETSc ordering.
99447c6ae99SBarry Smith   */
995fcfd50ebSBarry Smith   if (!da->prealloc_only) {
99647c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
99747c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
99847c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
99947c6ae99SBarry Smith 
10001321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10011321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
100247c6ae99SBarry Smith 
100347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
100447c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
100547c6ae99SBarry Smith 
10061321219cSEthan Coon         lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10071321219cSEthan Coon         lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
100847c6ae99SBarry Smith 
100947c6ae99SBarry Smith         for (k=0; k<nc; k++) {
101047c6ae99SBarry Smith           cnt = 0;
101147c6ae99SBarry Smith           for (l=lstart; l<lend+1; l++) {
101247c6ae99SBarry Smith             for (p=pstart; p<pend+1; p++) {
101347c6ae99SBarry Smith               if (l || p) {
1014aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
10158865f1eaSKarl Rupp                   for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p);
101647c6ae99SBarry Smith                 }
101747c6ae99SBarry Smith               } else {
101847c6ae99SBarry Smith                 if (dfill) {
10198865f1eaSKarl Rupp                   for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p);
102047c6ae99SBarry Smith                 } else {
10218865f1eaSKarl Rupp                   for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + gnx*l + p);
102247c6ae99SBarry Smith                 }
102347c6ae99SBarry Smith               }
102447c6ae99SBarry Smith             }
102547c6ae99SBarry Smith           }
102647c6ae99SBarry Smith           row  = k + nc*(slot);
102747c6ae99SBarry Smith           ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
102847c6ae99SBarry Smith         }
102947c6ae99SBarry Smith       }
103047c6ae99SBarry Smith     }
103147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
103247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
103347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
103447c6ae99SBarry Smith   }
103547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
103647c6ae99SBarry Smith   PetscFunctionReturn(0);
103747c6ae99SBarry Smith }
103847c6ae99SBarry Smith 
103947c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
104047c6ae99SBarry Smith 
104147c6ae99SBarry Smith #undef __FUNCT__
1042950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ"
1043950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM da,Mat J)
104447c6ae99SBarry Smith {
104547c6ae99SBarry Smith   PetscErrorCode         ierr;
104647c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
10470298fd71SBarry Smith   PetscInt               m,n,dim,s,*cols = NULL,k,nc,*rows = NULL,col,cnt,l,p,*dnz = NULL,*onz = NULL;
104847c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
104947c6ae99SBarry Smith   MPI_Comm               comm;
105047c6ae99SBarry Smith   PetscScalar            *values;
10511321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
105247c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1053aa219208SBarry Smith   DMDAStencilType        st;
105447c6ae99SBarry Smith 
105547c6ae99SBarry Smith   PetscFunctionBegin;
105647c6ae99SBarry Smith   /*
105747c6ae99SBarry Smith          nc - number of components per grid point
105847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
105947c6ae99SBarry Smith 
106047c6ae99SBarry Smith   */
10611321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
106247c6ae99SBarry Smith   col  = 2*s + 1;
106347c6ae99SBarry Smith 
1064aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1065aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
106647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
106747c6ae99SBarry Smith 
106847c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
10691411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
10701411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
107147c6ae99SBarry Smith 
107247c6ae99SBarry Smith   /* determine the matrix preallocation information */
107347c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
107447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
10751321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10761321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
107747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
10781321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10791321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
108047c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
10811321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
10821321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
108347c6ae99SBarry Smith 
108447c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
108547c6ae99SBarry Smith 
108647c6ae99SBarry Smith         cnt = 0;
108747c6ae99SBarry Smith         for (l=0; l<nc; l++) {
108847c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
108947c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
109047c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1091aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
109247c6ae99SBarry Smith                   cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
109347c6ae99SBarry Smith                 }
109447c6ae99SBarry Smith               }
109547c6ae99SBarry Smith             }
109647c6ae99SBarry Smith           }
109747c6ae99SBarry Smith           rows[l] = l + nc*(slot);
109847c6ae99SBarry Smith         }
1099784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
110047c6ae99SBarry Smith       }
110147c6ae99SBarry Smith     }
110247c6ae99SBarry Smith   }
1103f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
110447c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
110547c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
110647c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1107784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1108784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
110947c6ae99SBarry Smith 
111047c6ae99SBarry Smith   /*
111147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
111247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
111347c6ae99SBarry Smith     PETSc ordering.
111447c6ae99SBarry Smith   */
1115fcfd50ebSBarry Smith   if (!da->prealloc_only) {
111647c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
111747c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
111847c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
11191321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
11201321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
112147c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
11221321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
11231321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
112447c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
11251321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
11261321219cSEthan Coon           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
112747c6ae99SBarry Smith 
112847c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
112947c6ae99SBarry Smith 
113047c6ae99SBarry Smith           cnt = 0;
113147c6ae99SBarry Smith           for (l=0; l<nc; l++) {
113247c6ae99SBarry Smith             for (ii=istart; ii<iend+1; ii++) {
113347c6ae99SBarry Smith               for (jj=jstart; jj<jend+1; jj++) {
113447c6ae99SBarry Smith                 for (kk=kstart; kk<kend+1; kk++) {
1135aa219208SBarry Smith                   if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
113647c6ae99SBarry Smith                     cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
113747c6ae99SBarry Smith                   }
113847c6ae99SBarry Smith                 }
113947c6ae99SBarry Smith               }
114047c6ae99SBarry Smith             }
114147c6ae99SBarry Smith             rows[l] = l + nc*(slot);
114247c6ae99SBarry Smith           }
114347c6ae99SBarry Smith           ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
114447c6ae99SBarry Smith         }
114547c6ae99SBarry Smith       }
114647c6ae99SBarry Smith     }
114747c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
114847c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
114947c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
115047c6ae99SBarry Smith   }
115147c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
115247c6ae99SBarry Smith   PetscFunctionReturn(0);
115347c6ae99SBarry Smith }
115447c6ae99SBarry Smith 
115547c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
115647c6ae99SBarry Smith 
115747c6ae99SBarry Smith #undef __FUNCT__
1158ce308e1dSBarry Smith #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ_Fill"
1159ce308e1dSBarry Smith PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM da,Mat J)
1160ce308e1dSBarry Smith {
1161ce308e1dSBarry Smith   PetscErrorCode         ierr;
1162ce308e1dSBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1163ce308e1dSBarry Smith   PetscInt               xs,nx,i,j,gxs,gnx,row,k,l;
11640298fd71SBarry Smith   PetscInt               m,dim,s,*cols = NULL,nc,col,cnt,*ocols;
1165ce308e1dSBarry Smith   PetscInt               *ofill = dd->ofill;
1166ce308e1dSBarry Smith   PetscScalar            *values;
1167ce308e1dSBarry Smith   DMDABoundaryType       bx;
1168ce308e1dSBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1169ce308e1dSBarry Smith   PetscMPIInt            rank,size;
1170ce308e1dSBarry Smith 
1171ce308e1dSBarry Smith   PetscFunctionBegin;
1172*ce94432eSBarry Smith   if (dd->bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"With fill provided not implemented with periodic boundary conditions");
1173*ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
1174*ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr);
1175ce308e1dSBarry Smith 
1176ce308e1dSBarry Smith   /*
1177ce308e1dSBarry Smith          nc - number of components per grid point
1178ce308e1dSBarry Smith          col - number of colors needed in one direction for single component problem
1179ce308e1dSBarry Smith 
1180ce308e1dSBarry Smith   */
1181ce308e1dSBarry Smith   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
1182ce308e1dSBarry Smith   col  = 2*s + 1;
1183ce308e1dSBarry Smith 
1184ce308e1dSBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1185ce308e1dSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
1186ce308e1dSBarry Smith 
1187ce308e1dSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
1188ce308e1dSBarry Smith   ierr = PetscMalloc2(nx*nc,PetscInt,&cols,nx*nc,PetscInt,&ocols);CHKERRQ(ierr);
1189ce308e1dSBarry Smith   ierr = PetscMemzero(cols,nx*nc*sizeof(PetscInt));CHKERRQ(ierr);
1190ce308e1dSBarry Smith   ierr = PetscMemzero(ocols,nx*nc*sizeof(PetscInt));CHKERRQ(ierr);
1191ce308e1dSBarry Smith 
1192ce308e1dSBarry Smith   /*
1193ce308e1dSBarry Smith         note should be smaller for first and last process with no periodic
1194ce308e1dSBarry Smith         does not handle dfill
1195ce308e1dSBarry Smith   */
1196ce308e1dSBarry Smith   cnt = 0;
1197ce308e1dSBarry Smith   /* coupling with process to the left */
1198ce308e1dSBarry Smith   for (i=0; i<s; i++) {
1199ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
1200ce308e1dSBarry Smith       ocols[cnt] = ((!rank) ? 0 : (s - i)*(ofill[j+1] - ofill[j]));
1201ce308e1dSBarry Smith       cols[cnt]  = nc + (s + i)*(ofill[j+1] - ofill[j]);
1202ce308e1dSBarry Smith       cnt++;
1203ce308e1dSBarry Smith     }
1204ce308e1dSBarry Smith   }
1205ce308e1dSBarry Smith   for (i=s; i<nx-s; i++) {
1206ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
1207ce308e1dSBarry Smith       cols[cnt] = nc + 2*s*(ofill[j+1] - ofill[j]);
1208ce308e1dSBarry Smith       cnt++;
1209ce308e1dSBarry Smith     }
1210ce308e1dSBarry Smith   }
1211ce308e1dSBarry Smith   /* coupling with process to the right */
1212ce308e1dSBarry Smith   for (i=nx-s; i<nx; i++) {
1213ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
1214ce308e1dSBarry Smith       ocols[cnt] = ((rank == (size-1)) ? 0 : (i - nx + s + 1)*(ofill[j+1] - ofill[j]));
1215ce308e1dSBarry Smith       cols[cnt]  = nc + (s + nx - i - 1)*(ofill[j+1] - ofill[j]);
1216ce308e1dSBarry Smith       cnt++;
1217ce308e1dSBarry Smith     }
1218ce308e1dSBarry Smith   }
1219ce308e1dSBarry Smith 
1220ce308e1dSBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,cols);CHKERRQ(ierr);
1221ce308e1dSBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,cols,0,ocols);CHKERRQ(ierr);
1222ce308e1dSBarry Smith   ierr = PetscFree2(cols,ocols);CHKERRQ(ierr);
1223ce308e1dSBarry Smith 
1224ce308e1dSBarry Smith   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1225ce308e1dSBarry Smith   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1226ce308e1dSBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1227ce308e1dSBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
1228ce308e1dSBarry Smith 
1229ce308e1dSBarry Smith   /*
1230ce308e1dSBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1231ce308e1dSBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1232ce308e1dSBarry Smith     PETSc ordering.
1233ce308e1dSBarry Smith   */
1234ce308e1dSBarry Smith   if (!da->prealloc_only) {
1235ce308e1dSBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1236ce308e1dSBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1237ce308e1dSBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1238ce308e1dSBarry Smith 
1239ce308e1dSBarry Smith     row = xs*nc;
1240ce308e1dSBarry Smith     /* coupling with process to the left */
1241ce308e1dSBarry Smith     for (i=xs; i<xs+s; i++) {
1242ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1243ce308e1dSBarry Smith         cnt = 0;
1244ce308e1dSBarry Smith         if (rank) {
1245ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1246ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1247ce308e1dSBarry Smith           }
1248ce308e1dSBarry Smith         }
1249ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1250ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1251ce308e1dSBarry Smith         }
1252ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1253ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1254ce308e1dSBarry Smith         }
1255ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1256ce308e1dSBarry Smith         row++;
1257ce308e1dSBarry Smith       }
1258ce308e1dSBarry Smith     }
1259ce308e1dSBarry Smith     for (i=xs+s; i<xs+nx-s; i++) {
1260ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1261ce308e1dSBarry Smith         cnt = 0;
1262ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1263ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1264ce308e1dSBarry Smith         }
1265ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1266ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1267ce308e1dSBarry Smith         }
1268ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1269ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1270ce308e1dSBarry Smith         }
1271ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1272ce308e1dSBarry Smith         row++;
1273ce308e1dSBarry Smith       }
1274ce308e1dSBarry Smith     }
1275ce308e1dSBarry Smith     /* coupling with process to the right */
1276ce308e1dSBarry Smith     for (i=xs+nx-s; i<xs+nx; i++) {
1277ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1278ce308e1dSBarry Smith         cnt = 0;
1279ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1280ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1281ce308e1dSBarry Smith         }
1282ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1283ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1284ce308e1dSBarry Smith         }
1285ce308e1dSBarry Smith         if (rank < size-1) {
1286ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1287ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1288ce308e1dSBarry Smith           }
1289ce308e1dSBarry Smith         }
1290ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1291ce308e1dSBarry Smith         row++;
1292ce308e1dSBarry Smith       }
1293ce308e1dSBarry Smith     }
1294ce308e1dSBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1295ce308e1dSBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1296ce308e1dSBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1297ce308e1dSBarry Smith     ierr = PetscFree(cols);CHKERRQ(ierr);
1298ce308e1dSBarry Smith   }
1299ce308e1dSBarry Smith   PetscFunctionReturn(0);
1300ce308e1dSBarry Smith }
1301ce308e1dSBarry Smith 
1302ce308e1dSBarry Smith /* ---------------------------------------------------------------------------------*/
1303ce308e1dSBarry Smith 
1304ce308e1dSBarry Smith #undef __FUNCT__
1305950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ"
1306950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM da,Mat J)
130747c6ae99SBarry Smith {
130847c6ae99SBarry Smith   PetscErrorCode         ierr;
130947c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
13100298fd71SBarry Smith   PetscInt               m,dim,s,*cols = NULL,nc,*rows = NULL,col,cnt,l;
131147c6ae99SBarry Smith   PetscInt               istart,iend;
131247c6ae99SBarry Smith   PetscScalar            *values;
13131321219cSEthan Coon   DMDABoundaryType       bx;
131447c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
131547c6ae99SBarry Smith 
131647c6ae99SBarry Smith   PetscFunctionBegin;
131747c6ae99SBarry Smith   /*
131847c6ae99SBarry Smith          nc - number of components per grid point
131947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
132047c6ae99SBarry Smith 
132147c6ae99SBarry Smith   */
13221321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
132347c6ae99SBarry Smith   col  = 2*s + 1;
132447c6ae99SBarry Smith 
1325aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1326aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
132747c6ae99SBarry Smith 
1328f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
132947c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
133047c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
133147c6ae99SBarry Smith 
13321411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
13331411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1334784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1335784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
133647c6ae99SBarry Smith 
133747c6ae99SBarry Smith   /*
133847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
133947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
134047c6ae99SBarry Smith     PETSc ordering.
134147c6ae99SBarry Smith   */
1342fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1343ce308e1dSBarry Smith     ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
134447c6ae99SBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
134547c6ae99SBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
134647c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
134747c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
134847c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
134947c6ae99SBarry Smith       slot   = i - gxs;
135047c6ae99SBarry Smith 
135147c6ae99SBarry Smith       cnt = 0;
135247c6ae99SBarry Smith       for (l=0; l<nc; l++) {
135347c6ae99SBarry Smith         for (i1=istart; i1<iend+1; i1++) {
135447c6ae99SBarry Smith           cols[cnt++] = l + nc*(slot + i1);
135547c6ae99SBarry Smith         }
135647c6ae99SBarry Smith         rows[l] = l + nc*(slot);
135747c6ae99SBarry Smith       }
135847c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
135947c6ae99SBarry Smith     }
136047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
136147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
136247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
136347c6ae99SBarry Smith     ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
1364ce308e1dSBarry Smith   }
136547c6ae99SBarry Smith   PetscFunctionReturn(0);
136647c6ae99SBarry Smith }
136747c6ae99SBarry Smith 
136847c6ae99SBarry Smith #undef __FUNCT__
1369950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIBAIJ"
1370950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM da,Mat J)
137147c6ae99SBarry Smith {
137247c6ae99SBarry Smith   PetscErrorCode         ierr;
137347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
137447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
137547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
137647c6ae99SBarry Smith   MPI_Comm               comm;
137747c6ae99SBarry Smith   PetscScalar            *values;
13781321219cSEthan Coon   DMDABoundaryType       bx,by;
1379aa219208SBarry Smith   DMDAStencilType        st;
138047c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
138147c6ae99SBarry Smith 
138247c6ae99SBarry Smith   PetscFunctionBegin;
138347c6ae99SBarry Smith   /*
138447c6ae99SBarry Smith      nc - number of components per grid point
138547c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
138647c6ae99SBarry Smith   */
13871321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
138847c6ae99SBarry Smith   col  = 2*s + 1;
138947c6ae99SBarry Smith 
1390aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1391aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
139247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
139347c6ae99SBarry Smith 
139447c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
139547c6ae99SBarry Smith 
13961411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
13971411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
139847c6ae99SBarry Smith 
139947c6ae99SBarry Smith   /* determine the matrix preallocation information */
140047c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
140147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
14021321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14031321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
140447c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
14051321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14061321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
140747c6ae99SBarry Smith       slot   = i - gxs + gnx*(j - gys);
140847c6ae99SBarry Smith 
140947c6ae99SBarry Smith       /* Find block columns in block row */
141047c6ae99SBarry Smith       cnt = 0;
141147c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
141247c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1413aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
141447c6ae99SBarry Smith             cols[cnt++] = slot + ii + gnx*jj;
141547c6ae99SBarry Smith           }
141647c6ae99SBarry Smith         }
141747c6ae99SBarry Smith       }
1418784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
141947c6ae99SBarry Smith     }
142047c6ae99SBarry Smith   }
142147c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
142247c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
142347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
142447c6ae99SBarry Smith 
1425784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1426784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
142747c6ae99SBarry Smith 
142847c6ae99SBarry Smith   /*
142947c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
143047c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
143147c6ae99SBarry Smith     PETSc ordering.
143247c6ae99SBarry Smith   */
1433fcfd50ebSBarry Smith   if (!da->prealloc_only) {
143447c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
143547c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
143647c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
14371321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14381321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
143947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
14401321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14411321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
144247c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
144347c6ae99SBarry Smith         cnt  = 0;
144447c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
144547c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1446aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
144747c6ae99SBarry Smith               cols[cnt++] = slot + ii + gnx*jj;
144847c6ae99SBarry Smith             }
144947c6ae99SBarry Smith           }
145047c6ae99SBarry Smith         }
145147c6ae99SBarry Smith         ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
145247c6ae99SBarry Smith       }
145347c6ae99SBarry Smith     }
145447c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
145547c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
145647c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
145747c6ae99SBarry Smith   }
145847c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
145947c6ae99SBarry Smith   PetscFunctionReturn(0);
146047c6ae99SBarry Smith }
146147c6ae99SBarry Smith 
146247c6ae99SBarry Smith #undef __FUNCT__
1463950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIBAIJ"
1464950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM da,Mat J)
146547c6ae99SBarry Smith {
146647c6ae99SBarry Smith   PetscErrorCode         ierr;
146747c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
146847c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
146947c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
147047c6ae99SBarry Smith   MPI_Comm               comm;
147147c6ae99SBarry Smith   PetscScalar            *values;
14721321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1473aa219208SBarry Smith   DMDAStencilType        st;
147447c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
147547c6ae99SBarry Smith 
147647c6ae99SBarry Smith   PetscFunctionBegin;
147747c6ae99SBarry Smith   /*
147847c6ae99SBarry Smith          nc - number of components per grid point
147947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
148047c6ae99SBarry Smith 
148147c6ae99SBarry Smith   */
14821321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
148347c6ae99SBarry Smith   col  = 2*s + 1;
148447c6ae99SBarry Smith 
1485aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1486aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
148747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
148847c6ae99SBarry Smith 
148947c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
149047c6ae99SBarry Smith 
14911411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
14921411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
149347c6ae99SBarry Smith 
149447c6ae99SBarry Smith   /* determine the matrix preallocation information */
149547c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
149647c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
14971321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14981321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
149947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
15001321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15011321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
150247c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
15031321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15041321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
150547c6ae99SBarry Smith 
150647c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
150747c6ae99SBarry Smith 
150847c6ae99SBarry Smith         /* Find block columns in block row */
150947c6ae99SBarry Smith         cnt = 0;
151047c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
151147c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
151247c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1513aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
151447c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
151547c6ae99SBarry Smith               }
151647c6ae99SBarry Smith             }
151747c6ae99SBarry Smith           }
151847c6ae99SBarry Smith         }
1519784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
152047c6ae99SBarry Smith       }
152147c6ae99SBarry Smith     }
152247c6ae99SBarry Smith   }
152347c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
152447c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
152547c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
152647c6ae99SBarry Smith 
1527784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1528784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
152947c6ae99SBarry Smith 
153047c6ae99SBarry Smith   /*
153147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
153247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
153347c6ae99SBarry Smith     PETSc ordering.
153447c6ae99SBarry Smith   */
1535fcfd50ebSBarry Smith   if (!da->prealloc_only) {
153647c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
153747c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
153847c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
15391321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
15401321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
154147c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
15421321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15431321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
154447c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
15451321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15461321219cSEthan Coon           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
154747c6ae99SBarry Smith 
154847c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
154947c6ae99SBarry Smith 
155047c6ae99SBarry Smith           cnt = 0;
155147c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
155247c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
155347c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1554aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
155547c6ae99SBarry Smith                   cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
155647c6ae99SBarry Smith                 }
155747c6ae99SBarry Smith               }
155847c6ae99SBarry Smith             }
155947c6ae99SBarry Smith           }
156047c6ae99SBarry Smith           ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
156147c6ae99SBarry Smith         }
156247c6ae99SBarry Smith       }
156347c6ae99SBarry Smith     }
156447c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
156547c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156647c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
156747c6ae99SBarry Smith   }
156847c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
156947c6ae99SBarry Smith   PetscFunctionReturn(0);
157047c6ae99SBarry Smith }
157147c6ae99SBarry Smith 
157247c6ae99SBarry Smith #undef __FUNCT__
157347c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular"
157447c6ae99SBarry Smith /*
157547c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
157647c6ae99SBarry Smith   identify in the local ordering with periodic domain.
157747c6ae99SBarry Smith */
157847c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
157947c6ae99SBarry Smith {
158047c6ae99SBarry Smith   PetscErrorCode ierr;
158147c6ae99SBarry Smith   PetscInt       i,n;
158247c6ae99SBarry Smith 
158347c6ae99SBarry Smith   PetscFunctionBegin;
158447c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr);
158547c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr);
158647c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
158747c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
158847c6ae99SBarry Smith   }
158947c6ae99SBarry Smith   *cnt = n;
159047c6ae99SBarry Smith   PetscFunctionReturn(0);
159147c6ae99SBarry Smith }
159247c6ae99SBarry Smith 
159347c6ae99SBarry Smith #undef __FUNCT__
1594950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPISBAIJ"
1595950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM da,Mat J)
159647c6ae99SBarry Smith {
159747c6ae99SBarry Smith   PetscErrorCode         ierr;
159847c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
159947c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
160047c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
160147c6ae99SBarry Smith   MPI_Comm               comm;
160247c6ae99SBarry Smith   PetscScalar            *values;
16031321219cSEthan Coon   DMDABoundaryType       bx,by;
1604aa219208SBarry Smith   DMDAStencilType        st;
160547c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
160647c6ae99SBarry Smith 
160747c6ae99SBarry Smith   PetscFunctionBegin;
160847c6ae99SBarry Smith   /*
160947c6ae99SBarry Smith      nc - number of components per grid point
161047c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
161147c6ae99SBarry Smith   */
16121321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
161347c6ae99SBarry Smith   col  = 2*s + 1;
161447c6ae99SBarry Smith 
1615aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1616aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
161747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
161847c6ae99SBarry Smith 
161947c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
162047c6ae99SBarry Smith 
16211411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
16221411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
162347c6ae99SBarry Smith 
162447c6ae99SBarry Smith   /* determine the matrix preallocation information */
1625eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
162647c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
16271321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16281321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
162947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
16301321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16311321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
163247c6ae99SBarry Smith       slot   = i - gxs + gnx*(j - gys);
163347c6ae99SBarry Smith 
163447c6ae99SBarry Smith       /* Find block columns in block row */
163547c6ae99SBarry Smith       cnt = 0;
163647c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
163747c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1638aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) {
163947c6ae99SBarry Smith             cols[cnt++] = slot + ii + gnx*jj;
164047c6ae99SBarry Smith           }
164147c6ae99SBarry Smith         }
164247c6ae99SBarry Smith       }
164347c6ae99SBarry Smith       ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
164447c6ae99SBarry Smith       ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
164547c6ae99SBarry Smith     }
164647c6ae99SBarry Smith   }
164747c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
164847c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
164947c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
165047c6ae99SBarry Smith 
1651784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1652784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
165347c6ae99SBarry Smith 
165447c6ae99SBarry Smith   /*
165547c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
165647c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
165747c6ae99SBarry Smith     PETSc ordering.
165847c6ae99SBarry Smith   */
1659fcfd50ebSBarry Smith   if (!da->prealloc_only) {
166047c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
166147c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
166247c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
16631321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16641321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
166547c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
16661321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16671321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
166847c6ae99SBarry Smith         slot   = i - gxs + gnx*(j - gys);
166947c6ae99SBarry Smith 
167047c6ae99SBarry Smith         /* Find block columns in block row */
167147c6ae99SBarry Smith         cnt = 0;
167247c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
167347c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1674aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) {
167547c6ae99SBarry Smith               cols[cnt++] = slot + ii + gnx*jj;
167647c6ae99SBarry Smith             }
167747c6ae99SBarry Smith           }
167847c6ae99SBarry Smith         }
167947c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
168047c6ae99SBarry Smith         ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
168147c6ae99SBarry Smith       }
168247c6ae99SBarry Smith     }
168347c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
168447c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168547c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168647c6ae99SBarry Smith   }
168747c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
168847c6ae99SBarry Smith   PetscFunctionReturn(0);
168947c6ae99SBarry Smith }
169047c6ae99SBarry Smith 
169147c6ae99SBarry Smith #undef __FUNCT__
1692950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPISBAIJ"
1693950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM da,Mat J)
169447c6ae99SBarry Smith {
169547c6ae99SBarry Smith   PetscErrorCode         ierr;
169647c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
169747c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
169847c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
169947c6ae99SBarry Smith   MPI_Comm               comm;
170047c6ae99SBarry Smith   PetscScalar            *values;
17011321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1702aa219208SBarry Smith   DMDAStencilType        st;
170347c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
170447c6ae99SBarry Smith 
170547c6ae99SBarry Smith   PetscFunctionBegin;
170647c6ae99SBarry Smith   /*
170747c6ae99SBarry Smith      nc - number of components per grid point
170847c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
170947c6ae99SBarry Smith   */
17101321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
171147c6ae99SBarry Smith   col  = 2*s + 1;
171247c6ae99SBarry Smith 
1713aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1714aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
171547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
171647c6ae99SBarry Smith 
171747c6ae99SBarry Smith   /* create the matrix */
171847c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
171947c6ae99SBarry Smith 
17201411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
17211411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
172247c6ae99SBarry Smith 
172347c6ae99SBarry Smith   /* determine the matrix preallocation information */
1724eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
172547c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
17261321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
17271321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
172847c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
17291321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
17301321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
173147c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
17321321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
17331321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
173447c6ae99SBarry Smith 
173547c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
173647c6ae99SBarry Smith 
173747c6ae99SBarry Smith         /* Find block columns in block row */
173847c6ae99SBarry Smith         cnt = 0;
173947c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
174047c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
174147c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1742aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
174347c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
174447c6ae99SBarry Smith               }
174547c6ae99SBarry Smith             }
174647c6ae99SBarry Smith           }
174747c6ae99SBarry Smith         }
174847c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
174947c6ae99SBarry Smith         ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
175047c6ae99SBarry Smith       }
175147c6ae99SBarry Smith     }
175247c6ae99SBarry Smith   }
175347c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
175447c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
175547c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
175647c6ae99SBarry Smith 
1757784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1758784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
175947c6ae99SBarry Smith 
176047c6ae99SBarry Smith   /*
176147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
176247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
176347c6ae99SBarry Smith     PETSc ordering.
176447c6ae99SBarry Smith   */
1765fcfd50ebSBarry Smith   if (!da->prealloc_only) {
176647c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
176747c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
176847c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
17691321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
17701321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
177147c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
17721321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
17731321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
177447c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
17751321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
17761321219cSEthan Coon           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
177747c6ae99SBarry Smith 
177847c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
177947c6ae99SBarry Smith 
178047c6ae99SBarry Smith           cnt = 0;
178147c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
178247c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
178347c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1784aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
178547c6ae99SBarry Smith                   cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
178647c6ae99SBarry Smith                 }
178747c6ae99SBarry Smith               }
178847c6ae99SBarry Smith             }
178947c6ae99SBarry Smith           }
179047c6ae99SBarry Smith           ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
179147c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
179247c6ae99SBarry Smith         }
179347c6ae99SBarry Smith       }
179447c6ae99SBarry Smith     }
179547c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
179647c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
179747c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
179847c6ae99SBarry Smith   }
179947c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
180047c6ae99SBarry Smith   PetscFunctionReturn(0);
180147c6ae99SBarry Smith }
180247c6ae99SBarry Smith 
180347c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
180447c6ae99SBarry Smith 
180547c6ae99SBarry Smith #undef __FUNCT__
1806950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ_Fill"
1807950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J)
180847c6ae99SBarry Smith {
180947c6ae99SBarry Smith   PetscErrorCode         ierr;
181047c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
181147c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz;
181247c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
181347c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
181447c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
181547c6ae99SBarry Smith   MPI_Comm               comm;
181647c6ae99SBarry Smith   PetscScalar            *values;
18171321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
181847c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1819aa219208SBarry Smith   DMDAStencilType        st;
182047c6ae99SBarry Smith 
182147c6ae99SBarry Smith   PetscFunctionBegin;
182247c6ae99SBarry Smith   /*
182347c6ae99SBarry Smith          nc - number of components per grid point
182447c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
182547c6ae99SBarry Smith 
182647c6ae99SBarry Smith   */
18271321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
182847c6ae99SBarry Smith   col  = 2*s + 1;
1829*ce94432eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
183047c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1831*ce94432eSBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
183247c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1833*ce94432eSBarry Smith   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
183447c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
183547c6ae99SBarry Smith 
1836aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1837aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
183847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
183947c6ae99SBarry Smith 
184047c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
18411411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
18421411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
184347c6ae99SBarry Smith 
184447c6ae99SBarry Smith   /* determine the matrix preallocation information */
184547c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
184647c6ae99SBarry Smith 
184747c6ae99SBarry Smith 
184847c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
18491321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
18501321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
185147c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
18521321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
18531321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
185447c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
18551321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
18561321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
185747c6ae99SBarry Smith 
185847c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
185947c6ae99SBarry Smith 
186047c6ae99SBarry Smith         for (l=0; l<nc; l++) {
186147c6ae99SBarry Smith           cnt = 0;
186247c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
186347c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
186447c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
186547c6ae99SBarry Smith                 if (ii || jj || kk) {
1866aa219208SBarry Smith                   if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
18678865f1eaSKarl Rupp                     for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
186847c6ae99SBarry Smith                   }
186947c6ae99SBarry Smith                 } else {
187047c6ae99SBarry Smith                   if (dfill) {
18718865f1eaSKarl Rupp                     for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
187247c6ae99SBarry Smith                   } else {
18738865f1eaSKarl Rupp                     for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
187447c6ae99SBarry Smith                   }
187547c6ae99SBarry Smith                 }
187647c6ae99SBarry Smith               }
187747c6ae99SBarry Smith             }
187847c6ae99SBarry Smith           }
187947c6ae99SBarry Smith           row  = l + nc*(slot);
1880784ac674SJed Brown           ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
188147c6ae99SBarry Smith         }
188247c6ae99SBarry Smith       }
188347c6ae99SBarry Smith     }
188447c6ae99SBarry Smith   }
188547c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
188647c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
188747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1888784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1889784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
189047c6ae99SBarry Smith 
189147c6ae99SBarry Smith   /*
189247c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
189347c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
189447c6ae99SBarry Smith     PETSc ordering.
189547c6ae99SBarry Smith   */
1896fcfd50ebSBarry Smith   if (!da->prealloc_only) {
189747c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
189847c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
189947c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
19001321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
19011321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
190247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
19031321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
19041321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
190547c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
19061321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
19071321219cSEthan Coon           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
190847c6ae99SBarry Smith 
190947c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
191047c6ae99SBarry Smith 
191147c6ae99SBarry Smith           for (l=0; l<nc; l++) {
191247c6ae99SBarry Smith             cnt = 0;
191347c6ae99SBarry Smith             for (ii=istart; ii<iend+1; ii++) {
191447c6ae99SBarry Smith               for (jj=jstart; jj<jend+1; jj++) {
191547c6ae99SBarry Smith                 for (kk=kstart; kk<kend+1; kk++) {
191647c6ae99SBarry Smith                   if (ii || jj || kk) {
1917aa219208SBarry Smith                     if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
19188865f1eaSKarl Rupp                       for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
191947c6ae99SBarry Smith                     }
192047c6ae99SBarry Smith                   } else {
192147c6ae99SBarry Smith                     if (dfill) {
19228865f1eaSKarl Rupp                       for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
192347c6ae99SBarry Smith                     } else {
19248865f1eaSKarl Rupp                       for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
192547c6ae99SBarry Smith                     }
192647c6ae99SBarry Smith                   }
192747c6ae99SBarry Smith                 }
192847c6ae99SBarry Smith               }
192947c6ae99SBarry Smith             }
193047c6ae99SBarry Smith             row  = l + nc*(slot);
193147c6ae99SBarry Smith             ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
193247c6ae99SBarry Smith           }
193347c6ae99SBarry Smith         }
193447c6ae99SBarry Smith       }
193547c6ae99SBarry Smith     }
193647c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
193747c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
193847c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
193947c6ae99SBarry Smith   }
194047c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
194147c6ae99SBarry Smith   PetscFunctionReturn(0);
194247c6ae99SBarry Smith }
1943