xref: /petsc/src/dm/impls/da/fdda.c (revision b5579763ec3a855d2c812e379dcf34686eecb77d)
147c6ae99SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I      "petscdmda.h"     I*/
307475bc1SBarry Smith #include <petscmat.h>
447c6ae99SBarry Smith 
5e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring*);
6e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring*);
7e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring*);
8e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM,ISColoringType,ISColoring*);
947c6ae99SBarry Smith 
1047c6ae99SBarry Smith /*
1147c6ae99SBarry Smith    For ghost i that may be negative or greater than the upper bound this
1247c6ae99SBarry Smith   maps it into the 0:m-1 range using periodicity
1347c6ae99SBarry Smith */
1447c6ae99SBarry Smith #define SetInRange(i,m) ((i < 0) ? m+i : ((i >= m) ? i-m : i))
1547c6ae99SBarry Smith 
16ce308e1dSBarry Smith static PetscErrorCode DMDASetBlockFills_Private(const PetscInt *dfill,PetscInt w,PetscInt **rfill)
1747c6ae99SBarry Smith {
1847c6ae99SBarry Smith   PetscErrorCode ierr;
1947c6ae99SBarry Smith   PetscInt       i,j,nz,*fill;
2047c6ae99SBarry Smith 
2147c6ae99SBarry Smith   PetscFunctionBegin;
2247c6ae99SBarry Smith   if (!dfill) PetscFunctionReturn(0);
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   /* count number nonzeros */
2547c6ae99SBarry Smith   nz = 0;
2647c6ae99SBarry Smith   for (i=0; i<w; i++) {
2747c6ae99SBarry Smith     for (j=0; j<w; j++) {
2847c6ae99SBarry Smith       if (dfill[w*i+j]) nz++;
2947c6ae99SBarry Smith     }
3047c6ae99SBarry Smith   }
31854ce69bSBarry Smith   ierr = PetscMalloc1(nz + w + 1,&fill);CHKERRQ(ierr);
3247c6ae99SBarry Smith   /* construct modified CSR storage of nonzero structure */
33ce308e1dSBarry Smith   /*  fill[0 -- w] marks starts of each row of column indices (and end of last row)
34ce308e1dSBarry Smith    so fill[1] - fill[0] gives number of nonzeros in first row etc */
3547c6ae99SBarry Smith   nz = w + 1;
3647c6ae99SBarry Smith   for (i=0; i<w; i++) {
3747c6ae99SBarry Smith     fill[i] = nz;
3847c6ae99SBarry Smith     for (j=0; j<w; j++) {
3947c6ae99SBarry Smith       if (dfill[w*i+j]) {
4047c6ae99SBarry Smith         fill[nz] = j;
4147c6ae99SBarry Smith         nz++;
4247c6ae99SBarry Smith       }
4347c6ae99SBarry Smith     }
4447c6ae99SBarry Smith   }
4547c6ae99SBarry Smith   fill[w] = nz;
4647c6ae99SBarry Smith 
4747c6ae99SBarry Smith   *rfill = fill;
4847c6ae99SBarry Smith   PetscFunctionReturn(0);
4947c6ae99SBarry Smith }
5047c6ae99SBarry Smith 
5147c6ae99SBarry Smith /*@
52aa219208SBarry Smith     DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem
53950540a4SJed Brown     of the matrix returned by DMCreateMatrix().
5447c6ae99SBarry Smith 
55aa219208SBarry Smith     Logically Collective on DMDA
5647c6ae99SBarry Smith 
5747c6ae99SBarry Smith     Input Parameter:
5847c6ae99SBarry Smith +   da - the distributed array
590298fd71SBarry Smith .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
6047c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
6147c6ae99SBarry Smith 
6247c6ae99SBarry Smith 
6347c6ae99SBarry Smith     Level: developer
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith     Notes: This only makes sense when you are doing multicomponent problems but using the
6647c6ae99SBarry Smith        MPIAIJ matrix format
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
6947c6ae99SBarry Smith        representing coupling and 0 entries for missing coupling. For example
7047c6ae99SBarry Smith $             dfill[9] = {1, 0, 0,
7147c6ae99SBarry Smith $                         1, 1, 0,
7247c6ae99SBarry Smith $                         0, 1, 1}
7347c6ae99SBarry Smith        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
7447c6ae99SBarry Smith        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
7547c6ae99SBarry Smith        diagonal block).
7647c6ae99SBarry Smith 
77aa219208SBarry Smith      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
7847c6ae99SBarry Smith      can be represented in the dfill, ofill format
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith    Contributed by Glenn Hammond
8147c6ae99SBarry Smith 
828ddb5d8bSBarry Smith .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith @*/
85ce308e1dSBarry Smith PetscErrorCode  DMDASetBlockFills(DM da,const PetscInt *dfill,const PetscInt *ofill)
8647c6ae99SBarry Smith {
8747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
8847c6ae99SBarry Smith   PetscErrorCode ierr;
89ae4f298aSBarry Smith   PetscInt       i,k,cnt = 1;
9047c6ae99SBarry Smith 
9147c6ae99SBarry Smith   PetscFunctionBegin;
92aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr);
93aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr);
94ae4f298aSBarry Smith 
95ae4f298aSBarry Smith   /* ofillcount tracks the columns of ofill that have any nonzero in thems; the value in each location is the number of
96ae4f298aSBarry Smith    columns to the left with any nonzeros in them plus 1 */
971795a4d1SJed Brown   ierr = PetscCalloc1(dd->w,&dd->ofillcols);CHKERRQ(ierr);
98ae4f298aSBarry Smith   for (i=0; i<dd->w; i++) {
99ae4f298aSBarry Smith     for (k=dd->ofill[i]; k<dd->ofill[i+1]; k++) dd->ofillcols[dd->ofill[k]] = 1;
100ae4f298aSBarry Smith   }
101ae4f298aSBarry Smith   for (i=0; i<dd->w; i++) {
102ae4f298aSBarry Smith     if (dd->ofillcols[i]) {
103ae4f298aSBarry Smith       dd->ofillcols[i] = cnt++;
104ae4f298aSBarry Smith     }
105ae4f298aSBarry Smith   }
10647c6ae99SBarry Smith   PetscFunctionReturn(0);
10747c6ae99SBarry Smith }
10847c6ae99SBarry Smith 
10947c6ae99SBarry Smith 
110b412c318SBarry Smith PetscErrorCode  DMCreateColoring_DA(DM da,ISColoringType ctype,ISColoring *coloring)
11147c6ae99SBarry Smith {
11247c6ae99SBarry Smith   PetscErrorCode   ierr;
11347c6ae99SBarry Smith   PetscInt         dim,m,n,p,nc;
114bff4a2f0SMatthew G. Knepley   DMBoundaryType bx,by,bz;
11547c6ae99SBarry Smith   MPI_Comm         comm;
11647c6ae99SBarry Smith   PetscMPIInt      size;
11747c6ae99SBarry Smith   PetscBool        isBAIJ;
11847c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
11947c6ae99SBarry Smith 
12047c6ae99SBarry Smith   PetscFunctionBegin;
12147c6ae99SBarry Smith   /*
12247c6ae99SBarry Smith                                   m
12347c6ae99SBarry Smith           ------------------------------------------------------
12447c6ae99SBarry Smith          |                                                     |
12547c6ae99SBarry Smith          |                                                     |
12647c6ae99SBarry Smith          |               ----------------------                |
12747c6ae99SBarry Smith          |               |                    |                |
12847c6ae99SBarry Smith       n  |           yn  |                    |                |
12947c6ae99SBarry Smith          |               |                    |                |
13047c6ae99SBarry Smith          |               .---------------------                |
13147c6ae99SBarry Smith          |             (xs,ys)     xn                          |
13247c6ae99SBarry Smith          |            .                                        |
13347c6ae99SBarry Smith          |         (gxs,gys)                                   |
13447c6ae99SBarry Smith          |                                                     |
13547c6ae99SBarry Smith           -----------------------------------------------------
13647c6ae99SBarry Smith   */
13747c6ae99SBarry Smith 
13847c6ae99SBarry Smith   /*
13947c6ae99SBarry Smith          nc - number of components per grid point
14047c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
14147c6ae99SBarry Smith 
14247c6ae99SBarry Smith   */
1431321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&bx,&by,&bz,0);CHKERRQ(ierr);
14447c6ae99SBarry Smith 
14547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
14647c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1475bdb020cSBarry Smith   if (ctype == IS_COLORING_LOCAL) {
14847c6ae99SBarry Smith     if (size == 1) {
14947c6ae99SBarry Smith       ctype = IS_COLORING_GLOBAL;
15047c6ae99SBarry Smith     } else if (dim > 1) {
151bff4a2f0SMatthew G. Knepley       if ((m==1 && bx == DM_BOUNDARY_PERIODIC) || (n==1 && by == DM_BOUNDARY_PERIODIC) || (p==1 && bz == DM_BOUNDARY_PERIODIC)) {
1525bdb020cSBarry Smith         SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"IS_COLORING_LOCAL cannot be used for periodic boundary condition having both ends of the domain  on the same process");
15347c6ae99SBarry Smith       }
15447c6ae99SBarry Smith     }
15547c6ae99SBarry Smith   }
15647c6ae99SBarry Smith 
157aa219208SBarry Smith   /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ
15847c6ae99SBarry Smith      matrices is for the blocks, not the individual matrix elements  */
159b412c318SBarry Smith   ierr = PetscStrcmp(da->mattype,MATBAIJ,&isBAIJ);CHKERRQ(ierr);
160b412c318SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(da->mattype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);}
161b412c318SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(da->mattype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);}
16247c6ae99SBarry Smith   if (isBAIJ) {
16347c6ae99SBarry Smith     dd->w  = 1;
16447c6ae99SBarry Smith     dd->xs = dd->xs/nc;
16547c6ae99SBarry Smith     dd->xe = dd->xe/nc;
16647c6ae99SBarry Smith     dd->Xs = dd->Xs/nc;
16747c6ae99SBarry Smith     dd->Xe = dd->Xe/nc;
16847c6ae99SBarry Smith   }
16947c6ae99SBarry Smith 
17047c6ae99SBarry Smith   /*
171aa219208SBarry Smith      We do not provide a getcoloring function in the DMDA operations because
172aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
17347c6ae99SBarry Smith    more low-level then matrices.
17447c6ae99SBarry Smith   */
17547c6ae99SBarry Smith   if (dim == 1) {
176e727c939SJed Brown     ierr = DMCreateColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
17747c6ae99SBarry Smith   } else if (dim == 2) {
178e727c939SJed Brown     ierr =  DMCreateColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
17947c6ae99SBarry Smith   } else if (dim == 3) {
180e727c939SJed Brown     ierr =  DMCreateColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
181ce94432eSBarry 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);
18247c6ae99SBarry Smith   if (isBAIJ) {
18347c6ae99SBarry Smith     dd->w  = nc;
18447c6ae99SBarry Smith     dd->xs = dd->xs*nc;
18547c6ae99SBarry Smith     dd->xe = dd->xe*nc;
18647c6ae99SBarry Smith     dd->Xs = dd->Xs*nc;
18747c6ae99SBarry Smith     dd->Xe = dd->Xe*nc;
18847c6ae99SBarry Smith   }
18947c6ae99SBarry Smith   PetscFunctionReturn(0);
19047c6ae99SBarry Smith }
19147c6ae99SBarry Smith 
19247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
19347c6ae99SBarry Smith 
194e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
19547c6ae99SBarry Smith {
19647c6ae99SBarry Smith   PetscErrorCode   ierr;
19747c6ae99SBarry Smith   PetscInt         xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col;
19847c6ae99SBarry Smith   PetscInt         ncolors;
19947c6ae99SBarry Smith   MPI_Comm         comm;
200bff4a2f0SMatthew G. Knepley   DMBoundaryType bx,by;
201aa219208SBarry Smith   DMDAStencilType  st;
20247c6ae99SBarry Smith   ISColoringValue  *colors;
20347c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
20447c6ae99SBarry Smith 
20547c6ae99SBarry Smith   PetscFunctionBegin;
20647c6ae99SBarry Smith   /*
20747c6ae99SBarry Smith          nc - number of components per grid point
20847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
20947c6ae99SBarry Smith 
21047c6ae99SBarry Smith   */
2111321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
21247c6ae99SBarry Smith   col  = 2*s + 1;
213aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
214aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
21547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
21647c6ae99SBarry Smith 
21747c6ae99SBarry Smith   /* special case as taught to us by Paul Hovland */
218aa219208SBarry Smith   if (st == DMDA_STENCIL_STAR && s == 1) {
219e727c939SJed Brown     ierr = DMCreateColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
22047c6ae99SBarry Smith   } else {
22147c6ae99SBarry Smith 
222bff4a2f0SMatthew G. Knepley     if (bx == DM_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\
22347c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", m, col);
224bff4a2f0SMatthew G. Knepley     if (by == DM_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\
22547c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", n, col);
22647c6ae99SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
22747c6ae99SBarry Smith       if (!dd->localcoloring) {
228785e854fSJed Brown         ierr = PetscMalloc1(nc*nx*ny,&colors);CHKERRQ(ierr);
22947c6ae99SBarry Smith         ii   = 0;
23047c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
23147c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
23247c6ae99SBarry Smith             for (k=0; k<nc; k++) {
23347c6ae99SBarry Smith               colors[ii++] = k + nc*((i % col) + col*(j % col));
23447c6ae99SBarry Smith             }
23547c6ae99SBarry Smith           }
23647c6ae99SBarry Smith         }
23747c6ae99SBarry Smith         ncolors = nc + nc*(col-1 + col*(col-1));
238aaf3ff59SMatthew G. Knepley         ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,PETSC_OWN_POINTER,&dd->localcoloring);CHKERRQ(ierr);
23947c6ae99SBarry Smith       }
24047c6ae99SBarry Smith       *coloring = dd->localcoloring;
2415bdb020cSBarry Smith     } else if (ctype == IS_COLORING_LOCAL) {
24247c6ae99SBarry Smith       if (!dd->ghostedcoloring) {
243785e854fSJed Brown         ierr = PetscMalloc1(nc*gnx*gny,&colors);CHKERRQ(ierr);
24447c6ae99SBarry Smith         ii   = 0;
24547c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
24647c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
24747c6ae99SBarry Smith             for (k=0; k<nc; k++) {
24847c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
24947c6ae99SBarry Smith               colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
25047c6ae99SBarry Smith             }
25147c6ae99SBarry Smith           }
25247c6ae99SBarry Smith         }
25347c6ae99SBarry Smith         ncolors = nc + nc*(col - 1 + col*(col-1));
254aaf3ff59SMatthew G. Knepley         ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,PETSC_OWN_POINTER,&dd->ghostedcoloring);CHKERRQ(ierr);
25547c6ae99SBarry Smith         /* PetscIntView(ncolors,(PetscInt*)colors,0); */
25647c6ae99SBarry Smith 
2575bdb020cSBarry Smith         ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_LOCAL);CHKERRQ(ierr);
25847c6ae99SBarry Smith       }
25947c6ae99SBarry Smith       *coloring = dd->ghostedcoloring;
260ce94432eSBarry Smith     } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
26147c6ae99SBarry Smith   }
26247c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
26347c6ae99SBarry Smith   PetscFunctionReturn(0);
26447c6ae99SBarry Smith }
26547c6ae99SBarry Smith 
26647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
26747c6ae99SBarry Smith 
268e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
26947c6ae99SBarry Smith {
27047c6ae99SBarry Smith   PetscErrorCode   ierr;
27147c6ae99SBarry 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;
27247c6ae99SBarry Smith   PetscInt         ncolors;
27347c6ae99SBarry Smith   MPI_Comm         comm;
274bff4a2f0SMatthew G. Knepley   DMBoundaryType bx,by,bz;
275aa219208SBarry Smith   DMDAStencilType  st;
27647c6ae99SBarry Smith   ISColoringValue  *colors;
27747c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
27847c6ae99SBarry Smith 
27947c6ae99SBarry Smith   PetscFunctionBegin;
28047c6ae99SBarry Smith   /*
28147c6ae99SBarry Smith          nc - number of components per grid point
28247c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
28347c6ae99SBarry Smith 
28447c6ae99SBarry Smith   */
2851321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
28647c6ae99SBarry Smith   col  = 2*s + 1;
287bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
28847c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
289bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC && (n % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
29047c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
291bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC && (p % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
29247c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
29347c6ae99SBarry Smith 
294aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
295aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
29647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
29747c6ae99SBarry Smith 
29847c6ae99SBarry Smith   /* create the coloring */
29947c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
30047c6ae99SBarry Smith     if (!dd->localcoloring) {
301785e854fSJed Brown       ierr = PetscMalloc1(nc*nx*ny*nz,&colors);CHKERRQ(ierr);
30247c6ae99SBarry Smith       ii   = 0;
30347c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
30447c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
30547c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
30647c6ae99SBarry Smith             for (l=0; l<nc; l++) {
30747c6ae99SBarry Smith               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
30847c6ae99SBarry Smith             }
30947c6ae99SBarry Smith           }
31047c6ae99SBarry Smith         }
31147c6ae99SBarry Smith       }
31247c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
313aaf3ff59SMatthew G. Knepley       ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,PETSC_OWN_POINTER,&dd->localcoloring);CHKERRQ(ierr);
31447c6ae99SBarry Smith     }
31547c6ae99SBarry Smith     *coloring = dd->localcoloring;
3165bdb020cSBarry Smith   } else if (ctype == IS_COLORING_LOCAL) {
31747c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
318785e854fSJed Brown       ierr = PetscMalloc1(nc*gnx*gny*gnz,&colors);CHKERRQ(ierr);
31947c6ae99SBarry Smith       ii   = 0;
32047c6ae99SBarry Smith       for (k=gzs; k<gzs+gnz; k++) {
32147c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
32247c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
32347c6ae99SBarry Smith             for (l=0; l<nc; l++) {
32447c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
32547c6ae99SBarry Smith               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
32647c6ae99SBarry Smith             }
32747c6ae99SBarry Smith           }
32847c6ae99SBarry Smith         }
32947c6ae99SBarry Smith       }
33047c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
331aaf3ff59SMatthew G. Knepley       ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,PETSC_OWN_POINTER,&dd->ghostedcoloring);CHKERRQ(ierr);
3325bdb020cSBarry Smith       ierr    = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_LOCAL);CHKERRQ(ierr);
33347c6ae99SBarry Smith     }
33447c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
335ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
33647c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
33747c6ae99SBarry Smith   PetscFunctionReturn(0);
33847c6ae99SBarry Smith }
33947c6ae99SBarry Smith 
34047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
34147c6ae99SBarry Smith 
342e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
34347c6ae99SBarry Smith {
34447c6ae99SBarry Smith   PetscErrorCode   ierr;
34547c6ae99SBarry Smith   PetscInt         xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
34647c6ae99SBarry Smith   PetscInt         ncolors;
34747c6ae99SBarry Smith   MPI_Comm         comm;
348bff4a2f0SMatthew G. Knepley   DMBoundaryType bx;
34947c6ae99SBarry Smith   ISColoringValue  *colors;
35047c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
35147c6ae99SBarry Smith 
35247c6ae99SBarry Smith   PetscFunctionBegin;
35347c6ae99SBarry Smith   /*
35447c6ae99SBarry Smith          nc - number of components per grid point
35547c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
35647c6ae99SBarry Smith 
35747c6ae99SBarry Smith   */
3581321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
35947c6ae99SBarry Smith   col  = 2*s + 1;
36047c6ae99SBarry Smith 
361bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points %d is divisible\n\
36231e6f798SBarry Smith                                                           by 2*stencil_width + 1 %d\n",(int)m,(int)col);
36347c6ae99SBarry Smith 
364aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
365aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
36647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
36747c6ae99SBarry Smith 
36847c6ae99SBarry Smith   /* create the coloring */
36947c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
37047c6ae99SBarry Smith     if (!dd->localcoloring) {
371785e854fSJed Brown       ierr = PetscMalloc1(nc*nx,&colors);CHKERRQ(ierr);
372ae4f298aSBarry Smith       if (dd->ofillcols) {
373ae4f298aSBarry Smith         PetscInt tc = 0;
374ae4f298aSBarry Smith         for (i=0; i<nc; i++) tc += (PetscInt) (dd->ofillcols[i] > 0);
375ae4f298aSBarry Smith         i1 = 0;
376ae4f298aSBarry Smith         for (i=xs; i<xs+nx; i++) {
377ae4f298aSBarry Smith           for (l=0; l<nc; l++) {
378ae4f298aSBarry Smith             if (dd->ofillcols[l] && (i % col)) {
379ae4f298aSBarry Smith               colors[i1++] =  nc - 1 + tc*((i % col) - 1) + dd->ofillcols[l];
380ae4f298aSBarry Smith             } else {
381ae4f298aSBarry Smith               colors[i1++] = l;
382ae4f298aSBarry Smith             }
383ae4f298aSBarry Smith           }
384ae4f298aSBarry Smith         }
385ae4f298aSBarry Smith         ncolors = nc + 2*s*tc;
386ae4f298aSBarry Smith       } else {
38747c6ae99SBarry Smith         i1 = 0;
38847c6ae99SBarry Smith         for (i=xs; i<xs+nx; i++) {
38947c6ae99SBarry Smith           for (l=0; l<nc; l++) {
39047c6ae99SBarry Smith             colors[i1++] = l + nc*(i % col);
39147c6ae99SBarry Smith           }
39247c6ae99SBarry Smith         }
39347c6ae99SBarry Smith         ncolors = nc + nc*(col-1);
394ae4f298aSBarry Smith       }
395aaf3ff59SMatthew G. Knepley       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,PETSC_OWN_POINTER,&dd->localcoloring);CHKERRQ(ierr);
39647c6ae99SBarry Smith     }
39747c6ae99SBarry Smith     *coloring = dd->localcoloring;
3985bdb020cSBarry Smith   } else if (ctype == IS_COLORING_LOCAL) {
39947c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
400785e854fSJed Brown       ierr = PetscMalloc1(nc*gnx,&colors);CHKERRQ(ierr);
40147c6ae99SBarry Smith       i1   = 0;
40247c6ae99SBarry Smith       for (i=gxs; i<gxs+gnx; i++) {
40347c6ae99SBarry Smith         for (l=0; l<nc; l++) {
40447c6ae99SBarry Smith           /* the complicated stuff is to handle periodic boundaries */
40547c6ae99SBarry Smith           colors[i1++] = l + nc*(SetInRange(i,m) % col);
40647c6ae99SBarry Smith         }
40747c6ae99SBarry Smith       }
40847c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
409aaf3ff59SMatthew G. Knepley       ierr    = ISColoringCreate(comm,ncolors,nc*gnx,colors,PETSC_OWN_POINTER,&dd->ghostedcoloring);CHKERRQ(ierr);
4105bdb020cSBarry Smith       ierr    = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_LOCAL);CHKERRQ(ierr);
41147c6ae99SBarry Smith     }
41247c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
413ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
41447c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
41547c6ae99SBarry Smith   PetscFunctionReturn(0);
41647c6ae99SBarry Smith }
41747c6ae99SBarry Smith 
418e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
41947c6ae99SBarry Smith {
42047c6ae99SBarry Smith   PetscErrorCode   ierr;
42147c6ae99SBarry Smith   PetscInt         xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
42247c6ae99SBarry Smith   PetscInt         ncolors;
42347c6ae99SBarry Smith   MPI_Comm         comm;
424bff4a2f0SMatthew G. Knepley   DMBoundaryType bx,by;
42547c6ae99SBarry Smith   ISColoringValue  *colors;
42647c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
42747c6ae99SBarry Smith 
42847c6ae99SBarry Smith   PetscFunctionBegin;
42947c6ae99SBarry Smith   /*
43047c6ae99SBarry Smith          nc - number of components per grid point
43147c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
43247c6ae99SBarry Smith 
43347c6ae99SBarry Smith   */
4341321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);CHKERRQ(ierr);
435aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
436aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
43747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
43847c6ae99SBarry Smith 
439bff4a2f0SMatthew G. Knepley   if (bx == DM_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");
440bff4a2f0SMatthew G. Knepley   if (by == DM_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");
44147c6ae99SBarry Smith 
44247c6ae99SBarry Smith   /* create the coloring */
44347c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
44447c6ae99SBarry Smith     if (!dd->localcoloring) {
445785e854fSJed Brown       ierr = PetscMalloc1(nc*nx*ny,&colors);CHKERRQ(ierr);
44647c6ae99SBarry Smith       ii   = 0;
44747c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
44847c6ae99SBarry Smith         for (i=xs; i<xs+nx; i++) {
44947c6ae99SBarry Smith           for (k=0; k<nc; k++) {
45047c6ae99SBarry Smith             colors[ii++] = k + nc*((3*j+i) % 5);
45147c6ae99SBarry Smith           }
45247c6ae99SBarry Smith         }
45347c6ae99SBarry Smith       }
45447c6ae99SBarry Smith       ncolors = 5*nc;
455aaf3ff59SMatthew G. Knepley       ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,PETSC_OWN_POINTER,&dd->localcoloring);CHKERRQ(ierr);
45647c6ae99SBarry Smith     }
45747c6ae99SBarry Smith     *coloring = dd->localcoloring;
4585bdb020cSBarry Smith   } else if (ctype == IS_COLORING_LOCAL) {
45947c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
460785e854fSJed Brown       ierr = PetscMalloc1(nc*gnx*gny,&colors);CHKERRQ(ierr);
46147c6ae99SBarry Smith       ii = 0;
46247c6ae99SBarry Smith       for (j=gys; j<gys+gny; j++) {
46347c6ae99SBarry Smith         for (i=gxs; i<gxs+gnx; i++) {
46447c6ae99SBarry Smith           for (k=0; k<nc; k++) {
46547c6ae99SBarry Smith             colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
46647c6ae99SBarry Smith           }
46747c6ae99SBarry Smith         }
46847c6ae99SBarry Smith       }
46947c6ae99SBarry Smith       ncolors = 5*nc;
470aaf3ff59SMatthew G. Knepley       ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,PETSC_OWN_POINTER,&dd->ghostedcoloring);CHKERRQ(ierr);
4715bdb020cSBarry Smith       ierr    = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_LOCAL);CHKERRQ(ierr);
47247c6ae99SBarry Smith     }
47347c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
474ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
47547c6ae99SBarry Smith   PetscFunctionReturn(0);
47647c6ae99SBarry Smith }
47747c6ae99SBarry Smith 
47847c6ae99SBarry Smith /* =========================================================================== */
479950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM,Mat);
480ce308e1dSBarry Smith extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM,Mat);
481950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM,Mat);
482950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM,Mat);
483950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM,Mat);
484950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM,Mat);
485950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM,Mat);
486950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM,Mat);
487950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM,Mat);
488950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM,Mat);
48947c6ae99SBarry Smith 
4908bbdbebaSMatthew G Knepley /*@C
491c688c046SMatthew G Knepley    MatSetupDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith    Logically Collective on Mat
49447c6ae99SBarry Smith 
49547c6ae99SBarry Smith    Input Parameters:
49647c6ae99SBarry Smith +  mat - the matrix
49747c6ae99SBarry Smith -  da - the da
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith    Level: intermediate
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith @*/
502c688c046SMatthew G Knepley PetscErrorCode MatSetupDM(Mat mat,DM da)
50347c6ae99SBarry Smith {
50447c6ae99SBarry Smith   PetscErrorCode ierr;
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith   PetscFunctionBegin;
50747c6ae99SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
50847c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
509c688c046SMatthew G Knepley   ierr = PetscTryMethod(mat,"MatSetupDM_C",(Mat,DM),(mat,da));CHKERRQ(ierr);
51047c6ae99SBarry Smith   PetscFunctionReturn(0);
51147c6ae99SBarry Smith }
51247c6ae99SBarry Smith 
5137087cfbeSBarry Smith PetscErrorCode  MatView_MPI_DA(Mat A,PetscViewer viewer)
51447c6ae99SBarry Smith {
5159a42bb27SBarry Smith   DM                da;
51647c6ae99SBarry Smith   PetscErrorCode    ierr;
51747c6ae99SBarry Smith   const char        *prefix;
51847c6ae99SBarry Smith   Mat               Anatural;
51947c6ae99SBarry Smith   AO                ao;
52047c6ae99SBarry Smith   PetscInt          rstart,rend,*petsc,i;
52147c6ae99SBarry Smith   IS                is;
52247c6ae99SBarry Smith   MPI_Comm          comm;
52374388724SJed Brown   PetscViewerFormat format;
52447c6ae99SBarry Smith 
52547c6ae99SBarry Smith   PetscFunctionBegin;
52674388724SJed Brown   /* Check whether we are just printing info, in which case MatView() already viewed everything we wanted to view */
52774388724SJed Brown   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
52874388724SJed Brown   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
52974388724SJed Brown 
53047c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
531c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
532ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
53347c6ae99SBarry Smith 
534aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
53547c6ae99SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
536854ce69bSBarry Smith   ierr = PetscMalloc1(rend-rstart,&petsc);CHKERRQ(ierr);
53747c6ae99SBarry Smith   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
53847c6ae99SBarry Smith   ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr);
53947c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
54047c6ae99SBarry Smith 
54147c6ae99SBarry Smith   /* call viewer on natural ordering */
5427dae84e0SHong Zhang   ierr = MatCreateSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr);
543fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
54447c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr);
54547c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr);
54647c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr);
547539c167fSBarry Smith   ierr = MatView(Anatural,viewer);CHKERRQ(ierr);
548fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
54947c6ae99SBarry Smith   PetscFunctionReturn(0);
55047c6ae99SBarry Smith }
55147c6ae99SBarry Smith 
5527087cfbeSBarry Smith PetscErrorCode  MatLoad_MPI_DA(Mat A,PetscViewer viewer)
55347c6ae99SBarry Smith {
5549a42bb27SBarry Smith   DM             da;
55547c6ae99SBarry Smith   PetscErrorCode ierr;
55647c6ae99SBarry Smith   Mat            Anatural,Aapp;
55747c6ae99SBarry Smith   AO             ao;
558539c167fSBarry Smith   PetscInt       rstart,rend,*app,i,m,n,M,N;
55947c6ae99SBarry Smith   IS             is;
56047c6ae99SBarry Smith   MPI_Comm       comm;
56147c6ae99SBarry Smith 
56247c6ae99SBarry Smith   PetscFunctionBegin;
56347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
564c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
565ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
56647c6ae99SBarry Smith 
56747c6ae99SBarry Smith   /* Load the matrix in natural ordering */
568ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&Anatural);CHKERRQ(ierr);
56947c6ae99SBarry Smith   ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr);
570539c167fSBarry Smith   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
571539c167fSBarry Smith   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
572539c167fSBarry Smith   ierr = MatSetSizes(Anatural,m,n,M,N);CHKERRQ(ierr);
57347c6ae99SBarry Smith   ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr);
57447c6ae99SBarry Smith 
57547c6ae99SBarry Smith   /* Map natural ordering to application ordering and create IS */
576aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
57747c6ae99SBarry Smith   ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr);
578854ce69bSBarry Smith   ierr = PetscMalloc1(rend-rstart,&app);CHKERRQ(ierr);
57947c6ae99SBarry Smith   for (i=rstart; i<rend; i++) app[i-rstart] = i;
58047c6ae99SBarry Smith   ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr);
58147c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith   /* Do permutation and replace header */
5847dae84e0SHong Zhang   ierr = MatCreateSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr);
58528be2f97SBarry Smith   ierr = MatHeaderReplace(A,&Aapp);CHKERRQ(ierr);
586fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
587fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
58847c6ae99SBarry Smith   PetscFunctionReturn(0);
58947c6ae99SBarry Smith }
59047c6ae99SBarry Smith 
591b412c318SBarry Smith PetscErrorCode DMCreateMatrix_DA(DM da, Mat *J)
59247c6ae99SBarry Smith {
59347c6ae99SBarry Smith   PetscErrorCode ierr;
59447c6ae99SBarry Smith   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
59547c6ae99SBarry Smith   Mat            A;
59647c6ae99SBarry Smith   MPI_Comm       comm;
59719fd82e9SBarry Smith   MatType        Atype;
59837d0c07bSMatthew G Knepley   PetscSection   section, sectionGlobal;
5990298fd71SBarry Smith   void           (*aij)(void)=NULL,(*baij)(void)=NULL,(*sbaij)(void)=NULL;
600b412c318SBarry Smith   MatType        mtype;
60147c6ae99SBarry Smith   PetscMPIInt    size;
60247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
60347c6ae99SBarry Smith 
60447c6ae99SBarry Smith   PetscFunctionBegin;
605607a6623SBarry Smith   ierr = MatInitializePackage();CHKERRQ(ierr);
606b412c318SBarry Smith   mtype = da->mattype;
60747c6ae99SBarry Smith 
60837d0c07bSMatthew G Knepley   ierr = DMGetDefaultSection(da, &section);CHKERRQ(ierr);
60937d0c07bSMatthew G Knepley   if (section) {
61037d0c07bSMatthew G Knepley     PetscInt  bs = -1;
61137d0c07bSMatthew G Knepley     PetscInt  localSize;
61237d0c07bSMatthew G Knepley     PetscBool isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isSymmetric;
61337d0c07bSMatthew G Knepley 
61437d0c07bSMatthew G Knepley     ierr = DMGetDefaultGlobalSection(da, &sectionGlobal);CHKERRQ(ierr);
61537d0c07bSMatthew G Knepley     ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize);CHKERRQ(ierr);
616*b5579763SJed Brown     ierr = MatCreate(PetscObjectComm((PetscObject)da),&A);CHKERRQ(ierr);
617*b5579763SJed Brown     ierr = MatSetSizes(A,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
618*b5579763SJed Brown     ierr = MatSetType(A,mtype);CHKERRQ(ierr);
61937d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype,MATSHELL,&isShell);CHKERRQ(ierr);
62037d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype,MATBAIJ,&isBlock);CHKERRQ(ierr);
62137d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype,MATSEQBAIJ,&isSeqBlock);CHKERRQ(ierr);
62237d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype,MATMPIBAIJ,&isMPIBlock);CHKERRQ(ierr);
62337d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype,MATSBAIJ,&isSymBlock);CHKERRQ(ierr);
62437d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype,MATSEQSBAIJ,&isSymSeqBlock);CHKERRQ(ierr);
62537d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype,MATMPISBAIJ,&isSymMPIBlock);CHKERRQ(ierr);
62637d0c07bSMatthew G Knepley     /* Check for symmetric storage */
62737d0c07bSMatthew G Knepley     isSymmetric = (PetscBool) (isSymBlock || isSymSeqBlock || isSymMPIBlock);
62837d0c07bSMatthew G Knepley     if (isSymmetric) {
62937d0c07bSMatthew G Knepley       ierr = MatSetOption(*J, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);
63037d0c07bSMatthew G Knepley     }
63137d0c07bSMatthew G Knepley     if (!isShell) {
63237d0c07bSMatthew G Knepley       PetscInt *dnz, *onz, *dnzu, *onzu, bsLocal;
63337d0c07bSMatthew G Knepley 
63437d0c07bSMatthew G Knepley       if (bs < 0) {
63537d0c07bSMatthew G Knepley         if (isBlock || isSeqBlock || isMPIBlock || isSymBlock || isSymSeqBlock || isSymMPIBlock) {
63637d0c07bSMatthew G Knepley           PetscInt pStart, pEnd, p, dof;
63737d0c07bSMatthew G Knepley 
63837d0c07bSMatthew G Knepley           ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
63937d0c07bSMatthew G Knepley           for (p = pStart; p < pEnd; ++p) {
64037d0c07bSMatthew G Knepley             ierr = PetscSectionGetDof(sectionGlobal, p, &dof);CHKERRQ(ierr);
64137d0c07bSMatthew G Knepley             if (dof) {
64237d0c07bSMatthew G Knepley               bs = dof;
64337d0c07bSMatthew G Knepley               break;
64437d0c07bSMatthew G Knepley             }
64537d0c07bSMatthew G Knepley           }
64637d0c07bSMatthew G Knepley         } else {
64737d0c07bSMatthew G Knepley           bs = 1;
64837d0c07bSMatthew G Knepley         }
64937d0c07bSMatthew G Knepley         /* Must have same blocksize on all procs (some might have no points) */
65037d0c07bSMatthew G Knepley         bsLocal = bs;
651b2566f29SBarry Smith         ierr    = MPIU_Allreduce(&bsLocal, &bs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
65237d0c07bSMatthew G Knepley       }
6531795a4d1SJed Brown       ierr = PetscCalloc4(localSize/bs, &dnz, localSize/bs, &onz, localSize/bs, &dnzu, localSize/bs, &onzu);CHKERRQ(ierr);
654552f7358SJed Brown       /* ierr = DMPlexPreallocateOperator(dm, bs, section, sectionGlobal, dnz, onz, dnzu, onzu, *J, fillMatrix);CHKERRQ(ierr); */
65537d0c07bSMatthew G Knepley       ierr = PetscFree4(dnz, onz, dnzu, onzu);CHKERRQ(ierr);
65637d0c07bSMatthew G Knepley     }
65737d0c07bSMatthew G Knepley   }
65847c6ae99SBarry Smith   /*
65947c6ae99SBarry Smith                                   m
66047c6ae99SBarry Smith           ------------------------------------------------------
66147c6ae99SBarry Smith          |                                                     |
66247c6ae99SBarry Smith          |                                                     |
66347c6ae99SBarry Smith          |               ----------------------                |
66447c6ae99SBarry Smith          |               |                    |                |
66547c6ae99SBarry Smith       n  |           ny  |                    |                |
66647c6ae99SBarry Smith          |               |                    |                |
66747c6ae99SBarry Smith          |               .---------------------                |
66847c6ae99SBarry Smith          |             (xs,ys)     nx                          |
66947c6ae99SBarry Smith          |            .                                        |
67047c6ae99SBarry Smith          |         (gxs,gys)                                   |
67147c6ae99SBarry Smith          |                                                     |
67247c6ae99SBarry Smith           -----------------------------------------------------
67347c6ae99SBarry Smith   */
67447c6ae99SBarry Smith 
67547c6ae99SBarry Smith   /*
67647c6ae99SBarry Smith          nc - number of components per grid point
67747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
67847c6ae99SBarry Smith 
67947c6ae99SBarry Smith   */
680e30e807fSPeter Brune   M   = dd->M;
681e30e807fSPeter Brune   N   = dd->N;
682e30e807fSPeter Brune   P   = dd->P;
683c73cfb54SMatthew G. Knepley   dim = da->dim;
684e30e807fSPeter Brune   dof = dd->w;
685e30e807fSPeter Brune   /* ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); */
686aa219208SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
68747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
68847c6ae99SBarry Smith   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
68947c6ae99SBarry Smith   ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr);
690b412c318SBarry Smith   ierr = MatSetType(A,mtype);CHKERRQ(ierr);
69195ee5b0eSBarry Smith   ierr = MatSetDM(A,da);CHKERRQ(ierr);
692b06ff27eSHong Zhang   if (da->structure_only) {
693b06ff27eSHong Zhang     ierr = MatSetOption(A,MAT_STRUCTURE_ONLY,PETSC_TRUE);CHKERRQ(ierr);
694b06ff27eSHong Zhang   }
69547c6ae99SBarry Smith   ierr = MatGetType(A,&Atype);CHKERRQ(ierr);
69647c6ae99SBarry Smith   /*
697aa219208SBarry Smith      We do not provide a getmatrix function in the DMDA operations because
698aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
69947c6ae99SBarry Smith    more low-level than matrices. This is kind of cheating but, cause sometimes
700aa219208SBarry Smith    we think of DMDA has higher level than matrices.
70147c6ae99SBarry Smith 
70247c6ae99SBarry Smith      We could switch based on Atype (or mtype), but we do not since the
70347c6ae99SBarry Smith    specialized setting routines depend only the particular preallocation
70447c6ae99SBarry Smith    details of the matrix, not the type itself.
70547c6ae99SBarry Smith   */
70647c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
70747c6ae99SBarry Smith   if (!aij) {
70847c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
70947c6ae99SBarry Smith   }
71047c6ae99SBarry Smith   if (!aij) {
71147c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
71247c6ae99SBarry Smith     if (!baij) {
71347c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
71447c6ae99SBarry Smith     }
71547c6ae99SBarry Smith     if (!baij) {
71647c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
71747c6ae99SBarry Smith       if (!sbaij) {
71847c6ae99SBarry Smith         ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
71947c6ae99SBarry Smith       }
72047c6ae99SBarry Smith     }
72147c6ae99SBarry Smith   }
72247c6ae99SBarry Smith   if (aij) {
72347c6ae99SBarry Smith     if (dim == 1) {
724ce308e1dSBarry Smith       if (dd->ofill) {
725ce308e1dSBarry Smith         ierr = DMCreateMatrix_DA_1d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
726ce308e1dSBarry Smith       } else {
727950540a4SJed Brown         ierr = DMCreateMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr);
728ce308e1dSBarry Smith       }
72947c6ae99SBarry Smith     } else if (dim == 2) {
73047c6ae99SBarry Smith       if (dd->ofill) {
731950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
73247c6ae99SBarry Smith       } else {
733950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr);
73447c6ae99SBarry Smith       }
73547c6ae99SBarry Smith     } else if (dim == 3) {
73647c6ae99SBarry Smith       if (dd->ofill) {
737950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
73847c6ae99SBarry Smith       } else {
739950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr);
74047c6ae99SBarry Smith       }
74147c6ae99SBarry Smith     }
74247c6ae99SBarry Smith   } else if (baij) {
74347c6ae99SBarry Smith     if (dim == 2) {
744950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr);
74547c6ae99SBarry Smith     } else if (dim == 3) {
746950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr);
747ce94432eSBarry 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);
74847c6ae99SBarry Smith   } else if (sbaij) {
74947c6ae99SBarry Smith     if (dim == 2) {
750950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr);
75147c6ae99SBarry Smith     } else if (dim == 3) {
752950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr);
753ce94432eSBarry 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);
754869776cdSLisandro Dalcin   } else {
75545b6f7e9SBarry Smith     ISLocalToGlobalMapping ltog;
756869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
7572949035bSJed Brown     ierr = MatSetUp(A);CHKERRQ(ierr);
758869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMapping(A,ltog,ltog);CHKERRQ(ierr);
75947c6ae99SBarry Smith   }
760aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
76147c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
762c688c046SMatthew G Knepley   ierr = MatSetDM(A,da);CHKERRQ(ierr);
76347c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
76447c6ae99SBarry Smith   if (size > 1) {
76547c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
76647c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MPI_DA);CHKERRQ(ierr);
76747c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void))MatLoad_MPI_DA);CHKERRQ(ierr);
76847c6ae99SBarry Smith   }
769*b5579763SJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
77047c6ae99SBarry Smith   *J = A;
77147c6ae99SBarry Smith   PetscFunctionReturn(0);
77247c6ae99SBarry Smith }
77347c6ae99SBarry Smith 
77447c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
775950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM da,Mat J)
77647c6ae99SBarry Smith {
77747c6ae99SBarry Smith   PetscErrorCode         ierr;
778c1154cd5SBarry 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,M,N;
77947c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
78047c6ae99SBarry Smith   MPI_Comm               comm;
78147c6ae99SBarry Smith   PetscScalar            *values;
782bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
78345b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
784aa219208SBarry Smith   DMDAStencilType        st;
785c1154cd5SBarry Smith   PetscBool              removedups = PETSC_FALSE;
78647c6ae99SBarry Smith 
78747c6ae99SBarry Smith   PetscFunctionBegin;
78847c6ae99SBarry Smith   /*
78947c6ae99SBarry Smith          nc - number of components per grid point
79047c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith   */
793c1154cd5SBarry Smith   ierr = DMDAGetInfo(da,&dim,&m,&n,&M,&N,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
79447c6ae99SBarry Smith   col  = 2*s + 1;
795c1154cd5SBarry Smith   /*
796c1154cd5SBarry Smith        With one processor in periodic domains in a skinny dimension the code will label nonzero columns multiple times
797c1154cd5SBarry Smith        because of "wrapping" around the end of the domain hitting an entry already counted in the other direction.
798c1154cd5SBarry Smith   */
799c1154cd5SBarry Smith   if (M == 1 && 2*s >= m) removedups = PETSC_TRUE;
800c1154cd5SBarry Smith   if (N == 1 && 2*s >= n) removedups = PETSC_TRUE;
801aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
802aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
80347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
80447c6ae99SBarry Smith 
805dcca6d9dSJed Brown   ierr = PetscMalloc2(nc,&rows,col*col*nc*nc,&cols);CHKERRQ(ierr);
8061411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
80747c6ae99SBarry Smith 
80806ca8cadSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
80947c6ae99SBarry Smith   /* determine the matrix preallocation information */
81047c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
81147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
81247c6ae99SBarry Smith 
813bff4a2f0SMatthew G. Knepley     pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
814bff4a2f0SMatthew G. Knepley     pend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
81547c6ae99SBarry Smith 
81647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
81747c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
81847c6ae99SBarry Smith 
819bff4a2f0SMatthew G. Knepley       lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
820bff4a2f0SMatthew G. Knepley       lend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
82147c6ae99SBarry Smith 
82247c6ae99SBarry Smith       cnt = 0;
82347c6ae99SBarry Smith       for (k=0; k<nc; k++) {
82447c6ae99SBarry Smith         for (l=lstart; l<lend+1; l++) {
82547c6ae99SBarry Smith           for (p=pstart; p<pend+1; p++) {
826aa219208SBarry Smith             if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
82747c6ae99SBarry Smith               cols[cnt++] = k + nc*(slot + gnx*l + p);
82847c6ae99SBarry Smith             }
82947c6ae99SBarry Smith           }
83047c6ae99SBarry Smith         }
83147c6ae99SBarry Smith         rows[k] = k + nc*(slot);
83247c6ae99SBarry Smith       }
833c1154cd5SBarry Smith       if (removedups) {
834c1154cd5SBarry Smith         ierr = MatPreallocateSetLocalRemoveDups(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
835c1154cd5SBarry Smith       } else {
836784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
83747c6ae99SBarry Smith       }
83847c6ae99SBarry Smith     }
839c1154cd5SBarry Smith   }
840f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
84147c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
84247c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
84347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
84447c6ae99SBarry Smith 
845784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
84647c6ae99SBarry Smith 
84747c6ae99SBarry Smith   /*
84847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
84947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
85047c6ae99SBarry Smith     PETSc ordering.
85147c6ae99SBarry Smith   */
852fcfd50ebSBarry Smith   if (!da->prealloc_only) {
8531795a4d1SJed Brown     ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr);
85447c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
85547c6ae99SBarry Smith 
856bff4a2f0SMatthew G. Knepley       pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
857bff4a2f0SMatthew G. Knepley       pend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
85847c6ae99SBarry Smith 
85947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
86047c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
86147c6ae99SBarry Smith 
862bff4a2f0SMatthew G. Knepley         lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
863bff4a2f0SMatthew G. Knepley         lend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
86447c6ae99SBarry Smith 
86547c6ae99SBarry Smith         cnt = 0;
86647c6ae99SBarry Smith         for (k=0; k<nc; k++) {
86747c6ae99SBarry Smith           for (l=lstart; l<lend+1; l++) {
86847c6ae99SBarry Smith             for (p=pstart; p<pend+1; p++) {
869aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
87047c6ae99SBarry Smith                 cols[cnt++] = k + nc*(slot + gnx*l + p);
87147c6ae99SBarry Smith               }
87247c6ae99SBarry Smith             }
87347c6ae99SBarry Smith           }
87447c6ae99SBarry Smith           rows[k] = k + nc*(slot);
87547c6ae99SBarry Smith         }
87647c6ae99SBarry Smith         ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
87747c6ae99SBarry Smith       }
87847c6ae99SBarry Smith     }
87947c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
88047c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
88147c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
882189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
88347c6ae99SBarry Smith   }
88447c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
88547c6ae99SBarry Smith   PetscFunctionReturn(0);
88647c6ae99SBarry Smith }
88747c6ae99SBarry Smith 
888950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J)
88947c6ae99SBarry Smith {
89047c6ae99SBarry Smith   PetscErrorCode         ierr;
89147c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
892c1154cd5SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,maxcnt = 0,l,p,M,N;
89347c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
89447c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
89547c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
89647c6ae99SBarry Smith   MPI_Comm               comm;
89747c6ae99SBarry Smith   PetscScalar            *values;
898bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
89945b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
900aa219208SBarry Smith   DMDAStencilType        st;
901c1154cd5SBarry Smith   PetscBool              removedups = PETSC_FALSE;
90247c6ae99SBarry Smith 
90347c6ae99SBarry Smith   PetscFunctionBegin;
90447c6ae99SBarry Smith   /*
90547c6ae99SBarry Smith          nc - number of components per grid point
90647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
90747c6ae99SBarry Smith 
90847c6ae99SBarry Smith   */
909c1154cd5SBarry Smith   ierr = DMDAGetInfo(da,&dim,&m,&n,&M,&N,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
91047c6ae99SBarry Smith   col  = 2*s + 1;
911c1154cd5SBarry Smith   /*
912c1154cd5SBarry Smith        With one processor in periodic domains in a skinny dimension the code will label nonzero columns multiple times
913c1154cd5SBarry Smith        because of "wrapping" around the end of the domain hitting an entry already counted in the other direction.
914c1154cd5SBarry Smith   */
915c1154cd5SBarry Smith   if (M == 1 && 2*s >= m) removedups = PETSC_TRUE;
916c1154cd5SBarry Smith   if (N == 1 && 2*s >= n) removedups = PETSC_TRUE;
917aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
918aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
91947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
92047c6ae99SBarry Smith 
9214b26d1cfSBarry Smith   ierr = PetscMalloc1(col*col*nc,&cols);CHKERRQ(ierr);
9221411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
92347c6ae99SBarry Smith 
92406ca8cadSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
92547c6ae99SBarry Smith   /* determine the matrix preallocation information */
92647c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
92747c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
92847c6ae99SBarry Smith 
929bff4a2f0SMatthew G. Knepley     pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
930bff4a2f0SMatthew G. Knepley     pend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
93147c6ae99SBarry Smith 
93247c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
93347c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
93447c6ae99SBarry Smith 
935bff4a2f0SMatthew G. Knepley       lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
936bff4a2f0SMatthew G. Knepley       lend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
93747c6ae99SBarry Smith 
93847c6ae99SBarry Smith       for (k=0; k<nc; k++) {
93947c6ae99SBarry Smith         cnt = 0;
94047c6ae99SBarry Smith         for (l=lstart; l<lend+1; l++) {
94147c6ae99SBarry Smith           for (p=pstart; p<pend+1; p++) {
94247c6ae99SBarry Smith             if (l || p) {
943aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
9448865f1eaSKarl Rupp                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p);
94547c6ae99SBarry Smith               }
94647c6ae99SBarry Smith             } else {
94747c6ae99SBarry Smith               if (dfill) {
9488865f1eaSKarl Rupp                 for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p);
94947c6ae99SBarry Smith               } else {
9508865f1eaSKarl Rupp                 for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + gnx*l + p);
95147c6ae99SBarry Smith               }
95247c6ae99SBarry Smith             }
95347c6ae99SBarry Smith           }
95447c6ae99SBarry Smith         }
95547c6ae99SBarry Smith         row    = k + nc*(slot);
956c0ab637bSBarry Smith         maxcnt = PetscMax(maxcnt,cnt);
957c1154cd5SBarry Smith         if (removedups) {
958c1154cd5SBarry Smith           ierr   = MatPreallocateSetLocalRemoveDups(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
959c1154cd5SBarry Smith         } else {
960784ac674SJed Brown           ierr   = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
96147c6ae99SBarry Smith         }
96247c6ae99SBarry Smith       }
96347c6ae99SBarry Smith     }
964c1154cd5SBarry Smith   }
96547c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
96647c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
96747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
968784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
96947c6ae99SBarry Smith 
97047c6ae99SBarry Smith   /*
97147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
97247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
97347c6ae99SBarry Smith     PETSc ordering.
97447c6ae99SBarry Smith   */
975fcfd50ebSBarry Smith   if (!da->prealloc_only) {
976c0ab637bSBarry Smith     ierr = PetscCalloc1(maxcnt,&values);CHKERRQ(ierr);
97747c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
97847c6ae99SBarry Smith 
979bff4a2f0SMatthew G. Knepley       pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
980bff4a2f0SMatthew G. Knepley       pend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
98147c6ae99SBarry Smith 
98247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
98347c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
98447c6ae99SBarry Smith 
985bff4a2f0SMatthew G. Knepley         lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
986bff4a2f0SMatthew G. Knepley         lend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
98747c6ae99SBarry Smith 
98847c6ae99SBarry Smith         for (k=0; k<nc; k++) {
98947c6ae99SBarry Smith           cnt = 0;
99047c6ae99SBarry Smith           for (l=lstart; l<lend+1; l++) {
99147c6ae99SBarry Smith             for (p=pstart; p<pend+1; p++) {
99247c6ae99SBarry Smith               if (l || p) {
993aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
9948865f1eaSKarl Rupp                   for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p);
99547c6ae99SBarry Smith                 }
99647c6ae99SBarry Smith               } else {
99747c6ae99SBarry Smith                 if (dfill) {
9988865f1eaSKarl Rupp                   for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p);
99947c6ae99SBarry Smith                 } else {
10008865f1eaSKarl Rupp                   for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + gnx*l + p);
100147c6ae99SBarry Smith                 }
100247c6ae99SBarry Smith               }
100347c6ae99SBarry Smith             }
100447c6ae99SBarry Smith           }
100547c6ae99SBarry Smith           row  = k + nc*(slot);
100647c6ae99SBarry Smith           ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
100747c6ae99SBarry Smith         }
100847c6ae99SBarry Smith       }
100947c6ae99SBarry Smith     }
101047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
101147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1013189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
101447c6ae99SBarry Smith   }
101547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
101647c6ae99SBarry Smith   PetscFunctionReturn(0);
101747c6ae99SBarry Smith }
101847c6ae99SBarry Smith 
101947c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
102047c6ae99SBarry Smith 
1021950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM da,Mat J)
102247c6ae99SBarry Smith {
102347c6ae99SBarry Smith   PetscErrorCode         ierr;
102447c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
10250298fd71SBarry Smith   PetscInt               m,n,dim,s,*cols = NULL,k,nc,*rows = NULL,col,cnt,l,p,*dnz = NULL,*onz = NULL;
1026c1154cd5SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk,M,N,P;
102747c6ae99SBarry Smith   MPI_Comm               comm;
102847c6ae99SBarry Smith   PetscScalar            *values;
1029bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
103045b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
1031aa219208SBarry Smith   DMDAStencilType        st;
1032c1154cd5SBarry Smith   PetscBool              removedups = PETSC_FALSE;
103347c6ae99SBarry Smith 
103447c6ae99SBarry Smith   PetscFunctionBegin;
103547c6ae99SBarry Smith   /*
103647c6ae99SBarry Smith          nc - number of components per grid point
103747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
103847c6ae99SBarry Smith 
103947c6ae99SBarry Smith   */
1040c1154cd5SBarry Smith   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
104147c6ae99SBarry Smith   col  = 2*s + 1;
104247c6ae99SBarry Smith 
1043c1154cd5SBarry Smith   /*
1044c1154cd5SBarry Smith        With one processor in periodic domains in a skinny dimension the code will label nonzero columns multiple times
1045c1154cd5SBarry Smith        because of "wrapping" around the end of the domain hitting an entry already counted in the other direction.
1046c1154cd5SBarry Smith   */
1047c1154cd5SBarry Smith   if (M == 1 && 2*s >= m) removedups = PETSC_TRUE;
1048c1154cd5SBarry Smith   if (N == 1 && 2*s >= n) removedups = PETSC_TRUE;
1049c1154cd5SBarry Smith   if (P == 1 && 2*s >= p) removedups = PETSC_TRUE;
1050c1154cd5SBarry Smith 
1051aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1052aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
105347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
105447c6ae99SBarry Smith 
1055dcca6d9dSJed Brown   ierr = PetscMalloc2(nc,&rows,col*col*col*nc*nc,&cols);CHKERRQ(ierr);
10561411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
105747c6ae99SBarry Smith 
105806ca8cadSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
105947c6ae99SBarry Smith   /* determine the matrix preallocation information */
106047c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
106147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1062bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1063bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
106447c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1065bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1066bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
106747c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1068bff4a2f0SMatthew G. Knepley         kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1069bff4a2f0SMatthew G. Knepley         kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
107047c6ae99SBarry Smith 
107147c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
107247c6ae99SBarry Smith 
107347c6ae99SBarry Smith         cnt = 0;
107447c6ae99SBarry Smith         for (l=0; l<nc; l++) {
107547c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
107647c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
107747c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1078aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
107947c6ae99SBarry Smith                   cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
108047c6ae99SBarry Smith                 }
108147c6ae99SBarry Smith               }
108247c6ae99SBarry Smith             }
108347c6ae99SBarry Smith           }
108447c6ae99SBarry Smith           rows[l] = l + nc*(slot);
108547c6ae99SBarry Smith         }
1086c1154cd5SBarry Smith         if (removedups) {
1087c1154cd5SBarry Smith           ierr = MatPreallocateSetLocalRemoveDups(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
1088c1154cd5SBarry Smith         } else {
1089784ac674SJed Brown           ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
109047c6ae99SBarry Smith         }
109147c6ae99SBarry Smith       }
109247c6ae99SBarry Smith     }
1093c1154cd5SBarry Smith   }
1094f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
109547c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
109647c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
109747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1098784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
109947c6ae99SBarry Smith 
110047c6ae99SBarry Smith   /*
110147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
110247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
110347c6ae99SBarry Smith     PETSc ordering.
110447c6ae99SBarry Smith   */
1105fcfd50ebSBarry Smith   if (!da->prealloc_only) {
11061795a4d1SJed Brown     ierr = PetscCalloc1(col*col*col*nc*nc*nc,&values);CHKERRQ(ierr);
110747c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1108bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1109bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
111047c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1111bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1112bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
111347c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
1114bff4a2f0SMatthew G. Knepley           kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1115bff4a2f0SMatthew G. Knepley           kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
111647c6ae99SBarry Smith 
111747c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
111847c6ae99SBarry Smith 
111947c6ae99SBarry Smith           cnt = 0;
112047c6ae99SBarry Smith           for (l=0; l<nc; l++) {
112147c6ae99SBarry Smith             for (ii=istart; ii<iend+1; ii++) {
112247c6ae99SBarry Smith               for (jj=jstart; jj<jend+1; jj++) {
112347c6ae99SBarry Smith                 for (kk=kstart; kk<kend+1; kk++) {
1124aa219208SBarry Smith                   if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
112547c6ae99SBarry Smith                     cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
112647c6ae99SBarry Smith                   }
112747c6ae99SBarry Smith                 }
112847c6ae99SBarry Smith               }
112947c6ae99SBarry Smith             }
113047c6ae99SBarry Smith             rows[l] = l + nc*(slot);
113147c6ae99SBarry Smith           }
113247c6ae99SBarry Smith           ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
113347c6ae99SBarry Smith         }
113447c6ae99SBarry Smith       }
113547c6ae99SBarry Smith     }
113647c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
113747c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
113847c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1139189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
114047c6ae99SBarry Smith   }
114147c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
114247c6ae99SBarry Smith   PetscFunctionReturn(0);
114347c6ae99SBarry Smith }
114447c6ae99SBarry Smith 
114547c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
114647c6ae99SBarry Smith 
1147ce308e1dSBarry Smith PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM da,Mat J)
1148ce308e1dSBarry Smith {
1149ce308e1dSBarry Smith   PetscErrorCode         ierr;
1150ce308e1dSBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1151ce308e1dSBarry Smith   PetscInt               xs,nx,i,j,gxs,gnx,row,k,l;
11528d4c968fSBarry Smith   PetscInt               m,dim,s,*cols = NULL,nc,cnt,maxcnt = 0,*ocols;
11530acb5bebSBarry Smith   PetscInt               *ofill = dd->ofill,*dfill = dd->dfill;
1154ce308e1dSBarry Smith   PetscScalar            *values;
1155bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
115645b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
1157ce308e1dSBarry Smith   PetscMPIInt            rank,size;
1158ce308e1dSBarry Smith 
1159ce308e1dSBarry Smith   PetscFunctionBegin;
1160bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"With fill provided not implemented with periodic boundary conditions");
1161ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
1162ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr);
1163ce308e1dSBarry Smith 
1164ce308e1dSBarry Smith   /*
1165ce308e1dSBarry Smith          nc - number of components per grid point
1166ce308e1dSBarry Smith 
1167ce308e1dSBarry Smith   */
1168ce308e1dSBarry Smith   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
1169ce308e1dSBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1170ce308e1dSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
1171ce308e1dSBarry Smith 
1172ce308e1dSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
11731795a4d1SJed Brown   ierr = PetscCalloc2(nx*nc,&cols,nx*nc,&ocols);CHKERRQ(ierr);
1174ce308e1dSBarry Smith 
1175ce308e1dSBarry Smith   /*
1176ce308e1dSBarry Smith         note should be smaller for first and last process with no periodic
1177ce308e1dSBarry Smith         does not handle dfill
1178ce308e1dSBarry Smith   */
1179ce308e1dSBarry Smith   cnt = 0;
1180ce308e1dSBarry Smith   /* coupling with process to the left */
1181ce308e1dSBarry Smith   for (i=0; i<s; i++) {
1182ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
1183ce308e1dSBarry Smith       ocols[cnt] = ((!rank) ? 0 : (s - i)*(ofill[j+1] - ofill[j]));
11840acb5bebSBarry Smith       cols[cnt]  = dfill[j+1] - dfill[j] + (s + i)*(ofill[j+1] - ofill[j]);
1185c0ab637bSBarry Smith       maxcnt = PetscMax(maxcnt,ocols[cnt]+cols[cnt]);
1186ce308e1dSBarry Smith       cnt++;
1187ce308e1dSBarry Smith     }
1188ce308e1dSBarry Smith   }
1189ce308e1dSBarry Smith   for (i=s; i<nx-s; i++) {
1190ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
11910acb5bebSBarry Smith       cols[cnt] = dfill[j+1] - dfill[j] + 2*s*(ofill[j+1] - ofill[j]);
1192c0ab637bSBarry Smith       maxcnt = PetscMax(maxcnt,ocols[cnt]+cols[cnt]);
1193ce308e1dSBarry Smith       cnt++;
1194ce308e1dSBarry Smith     }
1195ce308e1dSBarry Smith   }
1196ce308e1dSBarry Smith   /* coupling with process to the right */
1197ce308e1dSBarry Smith   for (i=nx-s; i<nx; i++) {
1198ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
1199ce308e1dSBarry Smith       ocols[cnt] = ((rank == (size-1)) ? 0 : (i - nx + s + 1)*(ofill[j+1] - ofill[j]));
12000acb5bebSBarry Smith       cols[cnt]  = dfill[j+1] - dfill[j] + (s + nx - i - 1)*(ofill[j+1] - ofill[j]);
1201c0ab637bSBarry Smith       maxcnt = PetscMax(maxcnt,ocols[cnt]+cols[cnt]);
1202ce308e1dSBarry Smith       cnt++;
1203ce308e1dSBarry Smith     }
1204ce308e1dSBarry Smith   }
1205ce308e1dSBarry Smith 
1206ce308e1dSBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,cols);CHKERRQ(ierr);
1207ce308e1dSBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,cols,0,ocols);CHKERRQ(ierr);
1208ce308e1dSBarry Smith   ierr = PetscFree2(cols,ocols);CHKERRQ(ierr);
1209ce308e1dSBarry Smith 
1210ce308e1dSBarry Smith   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1211ce308e1dSBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1212ce308e1dSBarry Smith 
1213ce308e1dSBarry Smith   /*
1214ce308e1dSBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1215ce308e1dSBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1216ce308e1dSBarry Smith     PETSc ordering.
1217ce308e1dSBarry Smith   */
1218ce308e1dSBarry Smith   if (!da->prealloc_only) {
1219c0ab637bSBarry Smith     ierr = PetscCalloc2(maxcnt,&values,maxcnt,&cols);CHKERRQ(ierr);
1220ce308e1dSBarry Smith 
1221ce308e1dSBarry Smith     row = xs*nc;
1222ce308e1dSBarry Smith     /* coupling with process to the left */
1223ce308e1dSBarry Smith     for (i=xs; i<xs+s; i++) {
1224ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1225ce308e1dSBarry Smith         cnt = 0;
1226ce308e1dSBarry Smith         if (rank) {
1227ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1228ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1229ce308e1dSBarry Smith           }
1230ce308e1dSBarry Smith         }
12310acb5bebSBarry Smith         if (dfill) {
12320acb5bebSBarry Smith           for (k=dfill[j]; k<dfill[j+1]; k++) {
12330acb5bebSBarry Smith             cols[cnt++] = i*nc + dfill[k];
12340acb5bebSBarry Smith           }
12350acb5bebSBarry Smith         } else {
1236ce308e1dSBarry Smith           for (k=0; k<nc; k++) {
1237ce308e1dSBarry Smith             cols[cnt++] = i*nc + k;
1238ce308e1dSBarry Smith           }
12390acb5bebSBarry Smith         }
1240ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1241ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1242ce308e1dSBarry Smith         }
1243ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1244ce308e1dSBarry Smith         row++;
1245ce308e1dSBarry Smith       }
1246ce308e1dSBarry Smith     }
1247ce308e1dSBarry Smith     for (i=xs+s; i<xs+nx-s; i++) {
1248ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1249ce308e1dSBarry Smith         cnt = 0;
1250ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1251ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1252ce308e1dSBarry Smith         }
12530acb5bebSBarry Smith         if (dfill) {
12540acb5bebSBarry Smith           for (k=dfill[j]; k<dfill[j+1]; k++) {
12550acb5bebSBarry Smith             cols[cnt++] = i*nc + dfill[k];
12560acb5bebSBarry Smith           }
12570acb5bebSBarry Smith         } else {
1258ce308e1dSBarry Smith           for (k=0; k<nc; k++) {
1259ce308e1dSBarry Smith             cols[cnt++] = i*nc + k;
1260ce308e1dSBarry Smith           }
12610acb5bebSBarry Smith         }
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         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1266ce308e1dSBarry Smith         row++;
1267ce308e1dSBarry Smith       }
1268ce308e1dSBarry Smith     }
1269ce308e1dSBarry Smith     /* coupling with process to the right */
1270ce308e1dSBarry Smith     for (i=xs+nx-s; i<xs+nx; i++) {
1271ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1272ce308e1dSBarry Smith         cnt = 0;
1273ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1274ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1275ce308e1dSBarry Smith         }
12760acb5bebSBarry Smith         if (dfill) {
12770acb5bebSBarry Smith           for (k=dfill[j]; k<dfill[j+1]; k++) {
12780acb5bebSBarry Smith             cols[cnt++] = i*nc + dfill[k];
12790acb5bebSBarry Smith           }
12800acb5bebSBarry Smith         } else {
1281ce308e1dSBarry Smith           for (k=0; k<nc; k++) {
1282ce308e1dSBarry Smith             cols[cnt++] = i*nc + k;
1283ce308e1dSBarry Smith           }
12840acb5bebSBarry 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     }
1294c0ab637bSBarry Smith     ierr = PetscFree2(values,cols);CHKERRQ(ierr);
1295ce308e1dSBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1296ce308e1dSBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1297189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1298ce308e1dSBarry Smith   }
1299ce308e1dSBarry Smith   PetscFunctionReturn(0);
1300ce308e1dSBarry Smith }
1301ce308e1dSBarry Smith 
1302ce308e1dSBarry Smith /* ---------------------------------------------------------------------------------*/
1303ce308e1dSBarry Smith 
1304950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM da,Mat J)
130547c6ae99SBarry Smith {
130647c6ae99SBarry Smith   PetscErrorCode         ierr;
130747c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
13080298fd71SBarry Smith   PetscInt               m,dim,s,*cols = NULL,nc,*rows = NULL,col,cnt,l;
130947c6ae99SBarry Smith   PetscInt               istart,iend;
131047c6ae99SBarry Smith   PetscScalar            *values;
1311bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
131245b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
131347c6ae99SBarry Smith 
131447c6ae99SBarry Smith   PetscFunctionBegin;
131547c6ae99SBarry Smith   /*
131647c6ae99SBarry Smith          nc - number of components per grid point
131747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
131847c6ae99SBarry Smith 
131947c6ae99SBarry Smith   */
13201321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
132147c6ae99SBarry Smith   col  = 2*s + 1;
132247c6ae99SBarry Smith 
1323aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1324aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
132547c6ae99SBarry Smith 
1326f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
132747c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
132847c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
132947c6ae99SBarry Smith 
13301411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1331784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
133247c6ae99SBarry Smith 
133347c6ae99SBarry Smith   /*
133447c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
133547c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
133647c6ae99SBarry Smith     PETSc ordering.
133747c6ae99SBarry Smith   */
1338fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1339dcca6d9dSJed Brown     ierr = PetscMalloc2(nc,&rows,col*nc*nc,&cols);CHKERRQ(ierr);
13401795a4d1SJed Brown     ierr = PetscCalloc1(col*nc*nc,&values);CHKERRQ(ierr);
134147c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
134247c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
134347c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
134447c6ae99SBarry Smith       slot   = i - gxs;
134547c6ae99SBarry Smith 
134647c6ae99SBarry Smith       cnt = 0;
134747c6ae99SBarry Smith       for (l=0; l<nc; l++) {
134847c6ae99SBarry Smith         for (i1=istart; i1<iend+1; i1++) {
134947c6ae99SBarry Smith           cols[cnt++] = l + nc*(slot + i1);
135047c6ae99SBarry Smith         }
135147c6ae99SBarry Smith         rows[l] = l + nc*(slot);
135247c6ae99SBarry Smith       }
135347c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
135447c6ae99SBarry Smith     }
135547c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
135647c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
135747c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1358189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
135947c6ae99SBarry Smith     ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
1360ce308e1dSBarry Smith   }
136147c6ae99SBarry Smith   PetscFunctionReturn(0);
136247c6ae99SBarry Smith }
136347c6ae99SBarry Smith 
1364950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM da,Mat J)
136547c6ae99SBarry Smith {
136647c6ae99SBarry Smith   PetscErrorCode         ierr;
136747c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
136847c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
136947c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
137047c6ae99SBarry Smith   MPI_Comm               comm;
137147c6ae99SBarry Smith   PetscScalar            *values;
1372bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
1373aa219208SBarry Smith   DMDAStencilType        st;
137445b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
137547c6ae99SBarry Smith 
137647c6ae99SBarry Smith   PetscFunctionBegin;
137747c6ae99SBarry Smith   /*
137847c6ae99SBarry Smith      nc - number of components per grid point
137947c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
138047c6ae99SBarry Smith   */
13811321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
138247c6ae99SBarry Smith   col  = 2*s + 1;
138347c6ae99SBarry Smith 
1384aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1385aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
138647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
138747c6ae99SBarry Smith 
1388785e854fSJed Brown   ierr = PetscMalloc1(col*col*nc*nc,&cols);CHKERRQ(ierr);
138947c6ae99SBarry Smith 
13901411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
139147c6ae99SBarry Smith 
139247c6ae99SBarry Smith   /* determine the matrix preallocation information */
139347c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
139447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1395bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1396bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
139747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1398bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1399bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
140047c6ae99SBarry Smith       slot   = i - gxs + gnx*(j - gys);
140147c6ae99SBarry Smith 
140247c6ae99SBarry Smith       /* Find block columns in block row */
140347c6ae99SBarry Smith       cnt = 0;
140447c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
140547c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1406aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
140747c6ae99SBarry Smith             cols[cnt++] = slot + ii + gnx*jj;
140847c6ae99SBarry Smith           }
140947c6ae99SBarry Smith         }
141047c6ae99SBarry Smith       }
1411d6e23781SBarry Smith       ierr = MatPreallocateSetLocalBlock(ltog,1,&slot,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
141247c6ae99SBarry Smith     }
141347c6ae99SBarry Smith   }
141447c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
141547c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
141647c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
141747c6ae99SBarry Smith 
1418784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
141947c6ae99SBarry Smith 
142047c6ae99SBarry Smith   /*
142147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
142247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
142347c6ae99SBarry Smith     PETSc ordering.
142447c6ae99SBarry Smith   */
1425fcfd50ebSBarry Smith   if (!da->prealloc_only) {
14261795a4d1SJed Brown     ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr);
142747c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1428bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1429bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
143047c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1431bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1432bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
143347c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
143447c6ae99SBarry Smith         cnt  = 0;
143547c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
143647c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1437aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
143847c6ae99SBarry Smith               cols[cnt++] = slot + ii + gnx*jj;
143947c6ae99SBarry Smith             }
144047c6ae99SBarry Smith           }
144147c6ae99SBarry Smith         }
144247c6ae99SBarry Smith         ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
144347c6ae99SBarry Smith       }
144447c6ae99SBarry Smith     }
144547c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
144647c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
144747c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1448189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
144947c6ae99SBarry Smith   }
145047c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
145147c6ae99SBarry Smith   PetscFunctionReturn(0);
145247c6ae99SBarry Smith }
145347c6ae99SBarry Smith 
1454950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM da,Mat J)
145547c6ae99SBarry Smith {
145647c6ae99SBarry Smith   PetscErrorCode         ierr;
145747c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
145847c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
145947c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
146047c6ae99SBarry Smith   MPI_Comm               comm;
146147c6ae99SBarry Smith   PetscScalar            *values;
1462bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
1463aa219208SBarry Smith   DMDAStencilType        st;
146445b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
146547c6ae99SBarry Smith 
146647c6ae99SBarry Smith   PetscFunctionBegin;
146747c6ae99SBarry Smith   /*
146847c6ae99SBarry Smith          nc - number of components per grid point
146947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
147047c6ae99SBarry Smith 
147147c6ae99SBarry Smith   */
14721321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
147347c6ae99SBarry Smith   col  = 2*s + 1;
147447c6ae99SBarry Smith 
1475aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1476aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
147747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
147847c6ae99SBarry Smith 
1479785e854fSJed Brown   ierr = PetscMalloc1(col*col*col,&cols);CHKERRQ(ierr);
148047c6ae99SBarry Smith 
14811411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
148247c6ae99SBarry Smith 
148347c6ae99SBarry Smith   /* determine the matrix preallocation information */
148447c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
148547c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1486bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1487bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
148847c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1489bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1490bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
149147c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1492bff4a2f0SMatthew G. Knepley         kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1493bff4a2f0SMatthew G. Knepley         kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
149447c6ae99SBarry Smith 
149547c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
149647c6ae99SBarry Smith 
149747c6ae99SBarry Smith         /* Find block columns in block row */
149847c6ae99SBarry Smith         cnt = 0;
149947c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
150047c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
150147c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1502aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
150347c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
150447c6ae99SBarry Smith               }
150547c6ae99SBarry Smith             }
150647c6ae99SBarry Smith           }
150747c6ae99SBarry Smith         }
1508d6e23781SBarry Smith         ierr = MatPreallocateSetLocalBlock(ltog,1,&slot,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
150947c6ae99SBarry Smith       }
151047c6ae99SBarry Smith     }
151147c6ae99SBarry Smith   }
151247c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
151347c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
151447c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
151547c6ae99SBarry Smith 
1516784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
151747c6ae99SBarry Smith 
151847c6ae99SBarry Smith   /*
151947c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
152047c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
152147c6ae99SBarry Smith     PETSc ordering.
152247c6ae99SBarry Smith   */
1523fcfd50ebSBarry Smith   if (!da->prealloc_only) {
15241795a4d1SJed Brown     ierr = PetscCalloc1(col*col*col*nc*nc,&values);CHKERRQ(ierr);
152547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1526bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1527bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
152847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1529bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1530bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
153147c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
1532bff4a2f0SMatthew G. Knepley           kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1533bff4a2f0SMatthew G. Knepley           kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
153447c6ae99SBarry Smith 
153547c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
153647c6ae99SBarry Smith 
153747c6ae99SBarry Smith           cnt = 0;
153847c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
153947c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
154047c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1541aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
154247c6ae99SBarry Smith                   cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
154347c6ae99SBarry Smith                 }
154447c6ae99SBarry Smith               }
154547c6ae99SBarry Smith             }
154647c6ae99SBarry Smith           }
154747c6ae99SBarry Smith           ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
154847c6ae99SBarry Smith         }
154947c6ae99SBarry Smith       }
155047c6ae99SBarry Smith     }
155147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
155247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
155347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1554189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
155547c6ae99SBarry Smith   }
155647c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
155747c6ae99SBarry Smith   PetscFunctionReturn(0);
155847c6ae99SBarry Smith }
155947c6ae99SBarry Smith 
156047c6ae99SBarry Smith /*
156147c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
156247c6ae99SBarry Smith   identify in the local ordering with periodic domain.
156347c6ae99SBarry Smith */
156447c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
156547c6ae99SBarry Smith {
156647c6ae99SBarry Smith   PetscErrorCode ierr;
156747c6ae99SBarry Smith   PetscInt       i,n;
156847c6ae99SBarry Smith 
156947c6ae99SBarry Smith   PetscFunctionBegin;
1570d6e23781SBarry Smith   ierr = ISLocalToGlobalMappingApplyBlock(ltog,1,row,row);CHKERRQ(ierr);
1571d6e23781SBarry Smith   ierr = ISLocalToGlobalMappingApplyBlock(ltog,*cnt,col,col);CHKERRQ(ierr);
157247c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
157347c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
157447c6ae99SBarry Smith   }
157547c6ae99SBarry Smith   *cnt = n;
157647c6ae99SBarry Smith   PetscFunctionReturn(0);
157747c6ae99SBarry Smith }
157847c6ae99SBarry Smith 
1579950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM da,Mat J)
158047c6ae99SBarry Smith {
158147c6ae99SBarry Smith   PetscErrorCode         ierr;
158247c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
158347c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
158447c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
158547c6ae99SBarry Smith   MPI_Comm               comm;
158647c6ae99SBarry Smith   PetscScalar            *values;
1587bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
1588aa219208SBarry Smith   DMDAStencilType        st;
158945b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
159047c6ae99SBarry Smith 
159147c6ae99SBarry Smith   PetscFunctionBegin;
159247c6ae99SBarry Smith   /*
159347c6ae99SBarry Smith      nc - number of components per grid point
159447c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
159547c6ae99SBarry Smith   */
15961321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
159747c6ae99SBarry Smith   col  = 2*s + 1;
159847c6ae99SBarry Smith 
1599aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1600aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
160147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
160247c6ae99SBarry Smith 
1603785e854fSJed Brown   ierr = PetscMalloc1(col*col*nc*nc,&cols);CHKERRQ(ierr);
160447c6ae99SBarry Smith 
16051411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
160647c6ae99SBarry Smith 
160747c6ae99SBarry Smith   /* determine the matrix preallocation information */
1608eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
160947c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1610bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1611bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
161247c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1613bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1614bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
161547c6ae99SBarry Smith       slot   = i - gxs + gnx*(j - gys);
161647c6ae99SBarry Smith 
161747c6ae99SBarry Smith       /* Find block columns in block row */
161847c6ae99SBarry Smith       cnt = 0;
161947c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
162047c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1621aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) {
162247c6ae99SBarry Smith             cols[cnt++] = slot + ii + gnx*jj;
162347c6ae99SBarry Smith           }
162447c6ae99SBarry Smith         }
162547c6ae99SBarry Smith       }
162645b6f7e9SBarry Smith       ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr);
1627d6e23781SBarry Smith       ierr = MatPreallocateSymmetricSetBlock(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
162847c6ae99SBarry Smith     }
162947c6ae99SBarry Smith   }
163047c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
163147c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
163247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
163347c6ae99SBarry Smith 
1634784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
163547c6ae99SBarry Smith 
163647c6ae99SBarry Smith   /*
163747c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
163847c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
163947c6ae99SBarry Smith     PETSc ordering.
164047c6ae99SBarry Smith   */
1641fcfd50ebSBarry Smith   if (!da->prealloc_only) {
16421795a4d1SJed Brown     ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr);
164347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1644bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1645bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
164647c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1647bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1648bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
164947c6ae99SBarry Smith         slot   = i - gxs + gnx*(j - gys);
165047c6ae99SBarry Smith 
165147c6ae99SBarry Smith         /* Find block columns in block row */
165247c6ae99SBarry Smith         cnt = 0;
165347c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
165447c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1655aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) {
165647c6ae99SBarry Smith               cols[cnt++] = slot + ii + gnx*jj;
165747c6ae99SBarry Smith             }
165847c6ae99SBarry Smith           }
165947c6ae99SBarry Smith         }
166045b6f7e9SBarry Smith         ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr);
166147c6ae99SBarry Smith         ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
166247c6ae99SBarry Smith       }
166347c6ae99SBarry Smith     }
166447c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
166547c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166647c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1667189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
166847c6ae99SBarry Smith   }
166947c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
167047c6ae99SBarry Smith   PetscFunctionReturn(0);
167147c6ae99SBarry Smith }
167247c6ae99SBarry Smith 
1673950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM da,Mat J)
167447c6ae99SBarry Smith {
167547c6ae99SBarry Smith   PetscErrorCode         ierr;
167647c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
167747c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
167847c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
167947c6ae99SBarry Smith   MPI_Comm               comm;
168047c6ae99SBarry Smith   PetscScalar            *values;
1681bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
1682aa219208SBarry Smith   DMDAStencilType        st;
168345b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
168447c6ae99SBarry Smith 
168547c6ae99SBarry Smith   PetscFunctionBegin;
168647c6ae99SBarry Smith   /*
168747c6ae99SBarry Smith      nc - number of components per grid point
168847c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
168947c6ae99SBarry Smith   */
16901321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
169147c6ae99SBarry Smith   col  = 2*s + 1;
169247c6ae99SBarry Smith 
1693aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1694aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
169547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
169647c6ae99SBarry Smith 
169747c6ae99SBarry Smith   /* create the matrix */
1698785e854fSJed Brown   ierr = PetscMalloc1(col*col*col,&cols);CHKERRQ(ierr);
169947c6ae99SBarry Smith 
17001411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
170147c6ae99SBarry Smith 
170247c6ae99SBarry Smith   /* determine the matrix preallocation information */
1703eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
170447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1705bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1706bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
170747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1708bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1709bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
171047c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1711bff4a2f0SMatthew G. Knepley         kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1712bff4a2f0SMatthew G. Knepley         kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
171347c6ae99SBarry Smith 
171447c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
171547c6ae99SBarry Smith 
171647c6ae99SBarry Smith         /* Find block columns in block row */
171747c6ae99SBarry Smith         cnt = 0;
171847c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
171947c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
172047c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1721aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
172247c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
172347c6ae99SBarry Smith               }
172447c6ae99SBarry Smith             }
172547c6ae99SBarry Smith           }
172647c6ae99SBarry Smith         }
172745b6f7e9SBarry Smith         ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr);
1728d6e23781SBarry Smith         ierr = MatPreallocateSymmetricSetBlock(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
172947c6ae99SBarry Smith       }
173047c6ae99SBarry Smith     }
173147c6ae99SBarry Smith   }
173247c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
173347c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
173447c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
173547c6ae99SBarry Smith 
1736784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
173747c6ae99SBarry Smith 
173847c6ae99SBarry Smith   /*
173947c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
174047c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
174147c6ae99SBarry Smith     PETSc ordering.
174247c6ae99SBarry Smith   */
1743fcfd50ebSBarry Smith   if (!da->prealloc_only) {
17441795a4d1SJed Brown     ierr = PetscCalloc1(col*col*col*nc*nc,&values);CHKERRQ(ierr);
174547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1746bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1747bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
174847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1749bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1750bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
175147c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
1752bff4a2f0SMatthew G. Knepley           kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1753bff4a2f0SMatthew G. Knepley           kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
175447c6ae99SBarry Smith 
175547c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
175647c6ae99SBarry Smith 
175747c6ae99SBarry Smith           cnt = 0;
175847c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
175947c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
176047c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1761aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
176247c6ae99SBarry Smith                   cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
176347c6ae99SBarry Smith                 }
176447c6ae99SBarry Smith               }
176547c6ae99SBarry Smith             }
176647c6ae99SBarry Smith           }
176745b6f7e9SBarry Smith           ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr);
176847c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
176947c6ae99SBarry Smith         }
177047c6ae99SBarry Smith       }
177147c6ae99SBarry Smith     }
177247c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
177347c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177447c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1775189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
177647c6ae99SBarry Smith   }
177747c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
177847c6ae99SBarry Smith   PetscFunctionReturn(0);
177947c6ae99SBarry Smith }
178047c6ae99SBarry Smith 
178147c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
178247c6ae99SBarry Smith 
1783950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J)
178447c6ae99SBarry Smith {
178547c6ae99SBarry Smith   PetscErrorCode         ierr;
178647c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1787c0ab637bSBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt, maxcnt = 0,l,p,*dnz,*onz;
1788c1154cd5SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk,M,N,P;
178947c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
179047c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
179147c6ae99SBarry Smith   MPI_Comm               comm;
179247c6ae99SBarry Smith   PetscScalar            *values;
1793bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
179445b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
1795aa219208SBarry Smith   DMDAStencilType        st;
1796c1154cd5SBarry Smith   PetscBool              removedups = PETSC_FALSE;
179747c6ae99SBarry Smith 
179847c6ae99SBarry Smith   PetscFunctionBegin;
179947c6ae99SBarry Smith   /*
180047c6ae99SBarry Smith          nc - number of components per grid point
180147c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
180247c6ae99SBarry Smith 
180347c6ae99SBarry Smith   */
1804c1154cd5SBarry Smith   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
180547c6ae99SBarry Smith   col  = 2*s + 1;
1806bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
180747c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1808bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC && (n % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
180947c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1810bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC && (p % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
181147c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
181247c6ae99SBarry Smith 
1813c1154cd5SBarry Smith   /*
1814c1154cd5SBarry Smith        With one processor in periodic domains in a skinny dimension the code will label nonzero columns multiple times
1815c1154cd5SBarry Smith        because of "wrapping" around the end of the domain hitting an entry already counted in the other direction.
1816c1154cd5SBarry Smith   */
1817c1154cd5SBarry Smith   if (M == 1 && 2*s >= m) removedups = PETSC_TRUE;
1818c1154cd5SBarry Smith   if (N == 1 && 2*s >= n) removedups = PETSC_TRUE;
1819c1154cd5SBarry Smith   if (P == 1 && 2*s >= p) removedups = PETSC_TRUE;
1820c1154cd5SBarry Smith 
1821aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1822aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
182347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
182447c6ae99SBarry Smith 
1825785e854fSJed Brown   ierr = PetscMalloc1(col*col*col*nc,&cols);CHKERRQ(ierr);
18261411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
182747c6ae99SBarry Smith 
182847c6ae99SBarry Smith   /* determine the matrix preallocation information */
182947c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
183047c6ae99SBarry Smith 
183106ca8cadSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
183247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1833bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1834bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
183547c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1836bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1837bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
183847c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1839bff4a2f0SMatthew G. Knepley         kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1840bff4a2f0SMatthew G. Knepley         kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
184147c6ae99SBarry Smith 
184247c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
184347c6ae99SBarry Smith 
184447c6ae99SBarry Smith         for (l=0; l<nc; l++) {
184547c6ae99SBarry Smith           cnt = 0;
184647c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
184747c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
184847c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
184947c6ae99SBarry Smith                 if (ii || jj || kk) {
1850aa219208SBarry Smith                   if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
18518865f1eaSKarl 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);
185247c6ae99SBarry Smith                   }
185347c6ae99SBarry Smith                 } else {
185447c6ae99SBarry Smith                   if (dfill) {
18558865f1eaSKarl 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);
185647c6ae99SBarry Smith                   } else {
18578865f1eaSKarl Rupp                     for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
185847c6ae99SBarry Smith                   }
185947c6ae99SBarry Smith                 }
186047c6ae99SBarry Smith               }
186147c6ae99SBarry Smith             }
186247c6ae99SBarry Smith           }
186347c6ae99SBarry Smith           row  = l + nc*(slot);
1864c0ab637bSBarry Smith           maxcnt = PetscMax(maxcnt,cnt);
1865c1154cd5SBarry Smith           if (removedups) {
1866c1154cd5SBarry Smith             ierr = MatPreallocateSetLocalRemoveDups(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
1867c1154cd5SBarry Smith           } else {
1868784ac674SJed Brown             ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
186947c6ae99SBarry Smith           }
187047c6ae99SBarry Smith         }
187147c6ae99SBarry Smith       }
187247c6ae99SBarry Smith     }
1873c1154cd5SBarry Smith   }
187447c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
187547c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
187647c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1877784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
187847c6ae99SBarry Smith 
187947c6ae99SBarry Smith   /*
188047c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
188147c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
188247c6ae99SBarry Smith     PETSc ordering.
188347c6ae99SBarry Smith   */
1884fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1885c0ab637bSBarry Smith     ierr = PetscCalloc1(maxcnt,&values);CHKERRQ(ierr);
188647c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1887bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1888bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
188947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1890bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1891bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
189247c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
1893bff4a2f0SMatthew G. Knepley           kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1894bff4a2f0SMatthew G. Knepley           kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
189547c6ae99SBarry Smith 
189647c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
189747c6ae99SBarry Smith 
189847c6ae99SBarry Smith           for (l=0; l<nc; l++) {
189947c6ae99SBarry Smith             cnt = 0;
190047c6ae99SBarry Smith             for (ii=istart; ii<iend+1; ii++) {
190147c6ae99SBarry Smith               for (jj=jstart; jj<jend+1; jj++) {
190247c6ae99SBarry Smith                 for (kk=kstart; kk<kend+1; kk++) {
190347c6ae99SBarry Smith                   if (ii || jj || kk) {
1904aa219208SBarry Smith                     if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
19058865f1eaSKarl 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);
190647c6ae99SBarry Smith                     }
190747c6ae99SBarry Smith                   } else {
190847c6ae99SBarry Smith                     if (dfill) {
19098865f1eaSKarl 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);
191047c6ae99SBarry Smith                     } else {
19118865f1eaSKarl Rupp                       for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
191247c6ae99SBarry Smith                     }
191347c6ae99SBarry Smith                   }
191447c6ae99SBarry Smith                 }
191547c6ae99SBarry Smith               }
191647c6ae99SBarry Smith             }
191747c6ae99SBarry Smith             row  = l + nc*(slot);
191847c6ae99SBarry Smith             ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
191947c6ae99SBarry Smith           }
192047c6ae99SBarry Smith         }
192147c6ae99SBarry Smith       }
192247c6ae99SBarry Smith     }
192347c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
192447c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
192547c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1926189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
192747c6ae99SBarry Smith   }
192847c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
192947c6ae99SBarry Smith   PetscFunctionReturn(0);
193047c6ae99SBarry Smith }
1931