xref: /petsc/src/dm/impls/da/fdda.c (revision 1795a4d16c893ec2fc06bbbc6c5ce592a2de75d4)
147c6ae99SBarry Smith 
24035e84dSBarry Smith #include <petsc-private/dmdaimpl.h> /*I      "petscdmda.h"     I*/
307475bc1SBarry Smith #include <petscmat.h>
4b45d2f2cSJed Brown #include <petsc-private/matimpl.h>
547c6ae99SBarry Smith 
6e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring*);
7e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring*);
8e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring*);
9e727c939SJed Brown extern PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM,ISColoringType,ISColoring*);
1047c6ae99SBarry Smith 
1147c6ae99SBarry Smith /*
1247c6ae99SBarry Smith    For ghost i that may be negative or greater than the upper bound this
1347c6ae99SBarry Smith   maps it into the 0:m-1 range using periodicity
1447c6ae99SBarry Smith */
1547c6ae99SBarry Smith #define SetInRange(i,m) ((i < 0) ? m+i : ((i >= m) ? i-m : i))
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith #undef __FUNCT__
18aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills_Private"
19ce308e1dSBarry Smith static PetscErrorCode DMDASetBlockFills_Private(const PetscInt *dfill,PetscInt w,PetscInt **rfill)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   PetscErrorCode ierr;
2247c6ae99SBarry Smith   PetscInt       i,j,nz,*fill;
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   PetscFunctionBegin;
2547c6ae99SBarry Smith   if (!dfill) PetscFunctionReturn(0);
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   /* count number nonzeros */
2847c6ae99SBarry Smith   nz = 0;
2947c6ae99SBarry Smith   for (i=0; i<w; i++) {
3047c6ae99SBarry Smith     for (j=0; j<w; j++) {
3147c6ae99SBarry Smith       if (dfill[w*i+j]) nz++;
3247c6ae99SBarry Smith     }
3347c6ae99SBarry Smith   }
34785e854fSJed Brown   ierr = PetscMalloc1((nz + w + 1),&fill);CHKERRQ(ierr);
3547c6ae99SBarry Smith   /* construct modified CSR storage of nonzero structure */
36ce308e1dSBarry Smith   /*  fill[0 -- w] marks starts of each row of column indices (and end of last row)
37ce308e1dSBarry Smith    so fill[1] - fill[0] gives number of nonzeros in first row etc */
3847c6ae99SBarry Smith   nz = w + 1;
3947c6ae99SBarry Smith   for (i=0; i<w; i++) {
4047c6ae99SBarry Smith     fill[i] = nz;
4147c6ae99SBarry Smith     for (j=0; j<w; j++) {
4247c6ae99SBarry Smith       if (dfill[w*i+j]) {
4347c6ae99SBarry Smith         fill[nz] = j;
4447c6ae99SBarry Smith         nz++;
4547c6ae99SBarry Smith       }
4647c6ae99SBarry Smith     }
4747c6ae99SBarry Smith   }
4847c6ae99SBarry Smith   fill[w] = nz;
4947c6ae99SBarry Smith 
5047c6ae99SBarry Smith   *rfill = fill;
5147c6ae99SBarry Smith   PetscFunctionReturn(0);
5247c6ae99SBarry Smith }
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith #undef __FUNCT__
55aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills"
5647c6ae99SBarry Smith /*@
57aa219208SBarry Smith     DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem
58950540a4SJed Brown     of the matrix returned by DMCreateMatrix().
5947c6ae99SBarry Smith 
60aa219208SBarry Smith     Logically Collective on DMDA
6147c6ae99SBarry Smith 
6247c6ae99SBarry Smith     Input Parameter:
6347c6ae99SBarry Smith +   da - the distributed array
640298fd71SBarry Smith .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
6547c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
6647c6ae99SBarry Smith 
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith     Level: developer
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith     Notes: This only makes sense when you are doing multicomponent problems but using the
7147c6ae99SBarry Smith        MPIAIJ matrix format
7247c6ae99SBarry Smith 
7347c6ae99SBarry Smith            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
7447c6ae99SBarry Smith        representing coupling and 0 entries for missing coupling. For example
7547c6ae99SBarry Smith $             dfill[9] = {1, 0, 0,
7647c6ae99SBarry Smith $                         1, 1, 0,
7747c6ae99SBarry Smith $                         0, 1, 1}
7847c6ae99SBarry Smith        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
7947c6ae99SBarry Smith        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
8047c6ae99SBarry Smith        diagonal block).
8147c6ae99SBarry Smith 
82aa219208SBarry Smith      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
8347c6ae99SBarry Smith      can be represented in the dfill, ofill format
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith    Contributed by Glenn Hammond
8647c6ae99SBarry Smith 
878ddb5d8bSBarry Smith .seealso DMCreateMatrix(), DMDASetGetMatrix(), DMSetMatrixPreallocateOnly()
8847c6ae99SBarry Smith 
8947c6ae99SBarry Smith @*/
90ce308e1dSBarry Smith PetscErrorCode  DMDASetBlockFills(DM da,const PetscInt *dfill,const PetscInt *ofill)
9147c6ae99SBarry Smith {
9247c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
9347c6ae99SBarry Smith   PetscErrorCode ierr;
94ae4f298aSBarry Smith   PetscInt       i,k,cnt = 1;
9547c6ae99SBarry Smith 
9647c6ae99SBarry Smith   PetscFunctionBegin;
97aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr);
98aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr);
99ae4f298aSBarry Smith 
100ae4f298aSBarry Smith   /* ofillcount tracks the columns of ofill that have any nonzero in thems; the value in each location is the number of
101ae4f298aSBarry Smith    columns to the left with any nonzeros in them plus 1 */
102*1795a4d1SJed Brown   ierr = PetscCalloc1(dd->w,&dd->ofillcols);CHKERRQ(ierr);
103ae4f298aSBarry Smith   for (i=0; i<dd->w; i++) {
104ae4f298aSBarry Smith     for (k=dd->ofill[i]; k<dd->ofill[i+1]; k++) dd->ofillcols[dd->ofill[k]] = 1;
105ae4f298aSBarry Smith   }
106ae4f298aSBarry Smith   for (i=0; i<dd->w; i++) {
107ae4f298aSBarry Smith     if (dd->ofillcols[i]) {
108ae4f298aSBarry Smith       dd->ofillcols[i] = cnt++;
109ae4f298aSBarry Smith     }
110ae4f298aSBarry Smith   }
11147c6ae99SBarry Smith   PetscFunctionReturn(0);
11247c6ae99SBarry Smith }
11347c6ae99SBarry Smith 
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith #undef __FUNCT__
116e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA"
117b412c318SBarry Smith PetscErrorCode  DMCreateColoring_DA(DM da,ISColoringType ctype,ISColoring *coloring)
11847c6ae99SBarry Smith {
11947c6ae99SBarry Smith   PetscErrorCode   ierr;
12047c6ae99SBarry Smith   PetscInt         dim,m,n,p,nc;
1211321219cSEthan Coon   DMDABoundaryType bx,by,bz;
12247c6ae99SBarry Smith   MPI_Comm         comm;
12347c6ae99SBarry Smith   PetscMPIInt      size;
12447c6ae99SBarry Smith   PetscBool        isBAIJ;
12547c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith   PetscFunctionBegin;
12847c6ae99SBarry Smith   /*
12947c6ae99SBarry Smith                                   m
13047c6ae99SBarry Smith           ------------------------------------------------------
13147c6ae99SBarry Smith          |                                                     |
13247c6ae99SBarry Smith          |                                                     |
13347c6ae99SBarry Smith          |               ----------------------                |
13447c6ae99SBarry Smith          |               |                    |                |
13547c6ae99SBarry Smith       n  |           yn  |                    |                |
13647c6ae99SBarry Smith          |               |                    |                |
13747c6ae99SBarry Smith          |               .---------------------                |
13847c6ae99SBarry Smith          |             (xs,ys)     xn                          |
13947c6ae99SBarry Smith          |            .                                        |
14047c6ae99SBarry Smith          |         (gxs,gys)                                   |
14147c6ae99SBarry Smith          |                                                     |
14247c6ae99SBarry Smith           -----------------------------------------------------
14347c6ae99SBarry Smith   */
14447c6ae99SBarry Smith 
14547c6ae99SBarry Smith   /*
14647c6ae99SBarry Smith          nc - number of components per grid point
14747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
14847c6ae99SBarry Smith 
14947c6ae99SBarry Smith   */
1501321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&bx,&by,&bz,0);CHKERRQ(ierr);
15147c6ae99SBarry Smith 
15247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
15347c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
15447c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED) {
15547c6ae99SBarry Smith     if (size == 1) {
15647c6ae99SBarry Smith       ctype = IS_COLORING_GLOBAL;
15747c6ae99SBarry Smith     } else if (dim > 1) {
1581321219cSEthan Coon       if ((m==1 && bx == DMDA_BOUNDARY_PERIODIC) || (n==1 && by == DMDA_BOUNDARY_PERIODIC) || (p==1 && bz == DMDA_BOUNDARY_PERIODIC)) {
159ce94432eSBarry Smith         SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"IS_COLORING_GHOSTED cannot be used for periodic boundary condition having both ends of the domain  on the same process");
16047c6ae99SBarry Smith       }
16147c6ae99SBarry Smith     }
16247c6ae99SBarry Smith   }
16347c6ae99SBarry Smith 
164aa219208SBarry Smith   /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ
16547c6ae99SBarry Smith      matrices is for the blocks, not the individual matrix elements  */
166b412c318SBarry Smith   ierr = PetscStrcmp(da->mattype,MATBAIJ,&isBAIJ);CHKERRQ(ierr);
167b412c318SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(da->mattype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);}
168b412c318SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(da->mattype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);}
16947c6ae99SBarry Smith   if (isBAIJ) {
17047c6ae99SBarry Smith     dd->w  = 1;
17147c6ae99SBarry Smith     dd->xs = dd->xs/nc;
17247c6ae99SBarry Smith     dd->xe = dd->xe/nc;
17347c6ae99SBarry Smith     dd->Xs = dd->Xs/nc;
17447c6ae99SBarry Smith     dd->Xe = dd->Xe/nc;
17547c6ae99SBarry Smith   }
17647c6ae99SBarry Smith 
17747c6ae99SBarry Smith   /*
178aa219208SBarry Smith      We do not provide a getcoloring function in the DMDA operations because
179aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
18047c6ae99SBarry Smith    more low-level then matrices.
18147c6ae99SBarry Smith   */
18247c6ae99SBarry Smith   if (dim == 1) {
183e727c939SJed Brown     ierr = DMCreateColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
18447c6ae99SBarry Smith   } else if (dim == 2) {
185e727c939SJed Brown     ierr =  DMCreateColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
18647c6ae99SBarry Smith   } else if (dim == 3) {
187e727c939SJed Brown     ierr =  DMCreateColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
188ce94432eSBarry 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);
18947c6ae99SBarry Smith   if (isBAIJ) {
19047c6ae99SBarry Smith     dd->w  = nc;
19147c6ae99SBarry Smith     dd->xs = dd->xs*nc;
19247c6ae99SBarry Smith     dd->xe = dd->xe*nc;
19347c6ae99SBarry Smith     dd->Xs = dd->Xs*nc;
19447c6ae99SBarry Smith     dd->Xe = dd->Xe*nc;
19547c6ae99SBarry Smith   }
19647c6ae99SBarry Smith   PetscFunctionReturn(0);
19747c6ae99SBarry Smith }
19847c6ae99SBarry Smith 
19947c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
20047c6ae99SBarry Smith 
20147c6ae99SBarry Smith #undef __FUNCT__
202e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_MPIAIJ"
203e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
20447c6ae99SBarry Smith {
20547c6ae99SBarry Smith   PetscErrorCode   ierr;
20647c6ae99SBarry Smith   PetscInt         xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col;
20747c6ae99SBarry Smith   PetscInt         ncolors;
20847c6ae99SBarry Smith   MPI_Comm         comm;
2091321219cSEthan Coon   DMDABoundaryType bx,by;
210aa219208SBarry Smith   DMDAStencilType  st;
21147c6ae99SBarry Smith   ISColoringValue  *colors;
21247c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
21347c6ae99SBarry Smith 
21447c6ae99SBarry Smith   PetscFunctionBegin;
21547c6ae99SBarry Smith   /*
21647c6ae99SBarry Smith          nc - number of components per grid point
21747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
21847c6ae99SBarry Smith 
21947c6ae99SBarry Smith   */
2201321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
22147c6ae99SBarry Smith   col  = 2*s + 1;
222aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
223aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
22447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
22547c6ae99SBarry Smith 
22647c6ae99SBarry Smith   /* special case as taught to us by Paul Hovland */
227aa219208SBarry Smith   if (st == DMDA_STENCIL_STAR && s == 1) {
228e727c939SJed Brown     ierr = DMCreateColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
22947c6ae99SBarry Smith   } else {
23047c6ae99SBarry Smith 
231ce94432eSBarry Smith     if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\
23247c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", m, col);
233ce94432eSBarry Smith     if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\
23447c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", n, col);
23547c6ae99SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
23647c6ae99SBarry Smith       if (!dd->localcoloring) {
237785e854fSJed Brown         ierr = PetscMalloc1(nc*nx*ny,&colors);CHKERRQ(ierr);
23847c6ae99SBarry Smith         ii   = 0;
23947c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
24047c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
24147c6ae99SBarry Smith             for (k=0; k<nc; k++) {
24247c6ae99SBarry Smith               colors[ii++] = k + nc*((i % col) + col*(j % col));
24347c6ae99SBarry Smith             }
24447c6ae99SBarry Smith           }
24547c6ae99SBarry Smith         }
24647c6ae99SBarry Smith         ncolors = nc + nc*(col-1 + col*(col-1));
24747c6ae99SBarry Smith         ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
24847c6ae99SBarry Smith       }
24947c6ae99SBarry Smith       *coloring = dd->localcoloring;
25047c6ae99SBarry Smith     } else if (ctype == IS_COLORING_GHOSTED) {
25147c6ae99SBarry Smith       if (!dd->ghostedcoloring) {
252785e854fSJed Brown         ierr = PetscMalloc1(nc*gnx*gny,&colors);CHKERRQ(ierr);
25347c6ae99SBarry Smith         ii   = 0;
25447c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
25547c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
25647c6ae99SBarry Smith             for (k=0; k<nc; k++) {
25747c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
25847c6ae99SBarry Smith               colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
25947c6ae99SBarry Smith             }
26047c6ae99SBarry Smith           }
26147c6ae99SBarry Smith         }
26247c6ae99SBarry Smith         ncolors = nc + nc*(col - 1 + col*(col-1));
26347c6ae99SBarry Smith         ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
26447c6ae99SBarry Smith         /* PetscIntView(ncolors,(PetscInt*)colors,0); */
26547c6ae99SBarry Smith 
26647c6ae99SBarry Smith         ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
26747c6ae99SBarry Smith       }
26847c6ae99SBarry Smith       *coloring = dd->ghostedcoloring;
269ce94432eSBarry Smith     } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
27047c6ae99SBarry Smith   }
27147c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
27247c6ae99SBarry Smith   PetscFunctionReturn(0);
27347c6ae99SBarry Smith }
27447c6ae99SBarry Smith 
27547c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
27647c6ae99SBarry Smith 
27747c6ae99SBarry Smith #undef __FUNCT__
278e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_3d_MPIAIJ"
279e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
28047c6ae99SBarry Smith {
28147c6ae99SBarry Smith   PetscErrorCode   ierr;
28247c6ae99SBarry 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;
28347c6ae99SBarry Smith   PetscInt         ncolors;
28447c6ae99SBarry Smith   MPI_Comm         comm;
2851321219cSEthan Coon   DMDABoundaryType bx,by,bz;
286aa219208SBarry Smith   DMDAStencilType  st;
28747c6ae99SBarry Smith   ISColoringValue  *colors;
28847c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
28947c6ae99SBarry Smith 
29047c6ae99SBarry Smith   PetscFunctionBegin;
29147c6ae99SBarry Smith   /*
29247c6ae99SBarry Smith          nc - number of components per grid point
29347c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
29447c6ae99SBarry Smith 
29547c6ae99SBarry Smith   */
2961321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
29747c6ae99SBarry Smith   col  = 2*s + 1;
298ce94432eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
29947c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
300ce94432eSBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
30147c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
302ce94432eSBarry Smith   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
30347c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
30447c6ae99SBarry Smith 
305aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
306aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
30747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
30847c6ae99SBarry Smith 
30947c6ae99SBarry Smith   /* create the coloring */
31047c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
31147c6ae99SBarry Smith     if (!dd->localcoloring) {
312785e854fSJed Brown       ierr = PetscMalloc1(nc*nx*ny*nz,&colors);CHKERRQ(ierr);
31347c6ae99SBarry Smith       ii   = 0;
31447c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
31547c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
31647c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
31747c6ae99SBarry Smith             for (l=0; l<nc; l++) {
31847c6ae99SBarry Smith               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
31947c6ae99SBarry Smith             }
32047c6ae99SBarry Smith           }
32147c6ae99SBarry Smith         }
32247c6ae99SBarry Smith       }
32347c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
32447c6ae99SBarry Smith       ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr);
32547c6ae99SBarry Smith     }
32647c6ae99SBarry Smith     *coloring = dd->localcoloring;
32747c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
32847c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
329785e854fSJed Brown       ierr = PetscMalloc1(nc*gnx*gny*gnz,&colors);CHKERRQ(ierr);
33047c6ae99SBarry Smith       ii   = 0;
33147c6ae99SBarry Smith       for (k=gzs; k<gzs+gnz; k++) {
33247c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
33347c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
33447c6ae99SBarry Smith             for (l=0; l<nc; l++) {
33547c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
33647c6ae99SBarry Smith               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
33747c6ae99SBarry Smith             }
33847c6ae99SBarry Smith           }
33947c6ae99SBarry Smith         }
34047c6ae99SBarry Smith       }
34147c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
34247c6ae99SBarry Smith       ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
34347c6ae99SBarry Smith       ierr    = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
34447c6ae99SBarry Smith     }
34547c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
346ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
34747c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
34847c6ae99SBarry Smith   PetscFunctionReturn(0);
34947c6ae99SBarry Smith }
35047c6ae99SBarry Smith 
35147c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
35247c6ae99SBarry Smith 
35347c6ae99SBarry Smith #undef __FUNCT__
354e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_1d_MPIAIJ"
355e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
35647c6ae99SBarry Smith {
35747c6ae99SBarry Smith   PetscErrorCode   ierr;
35847c6ae99SBarry Smith   PetscInt         xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
35947c6ae99SBarry Smith   PetscInt         ncolors;
36047c6ae99SBarry Smith   MPI_Comm         comm;
3611321219cSEthan Coon   DMDABoundaryType bx;
36247c6ae99SBarry Smith   ISColoringValue  *colors;
36347c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
36447c6ae99SBarry Smith 
36547c6ae99SBarry Smith   PetscFunctionBegin;
36647c6ae99SBarry Smith   /*
36747c6ae99SBarry Smith          nc - number of components per grid point
36847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
36947c6ae99SBarry Smith 
37047c6ae99SBarry Smith   */
3711321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
37247c6ae99SBarry Smith   col  = 2*s + 1;
37347c6ae99SBarry Smith 
374ce94432eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points %d is divisible\n\
37531e6f798SBarry Smith                                                           by 2*stencil_width + 1 %d\n",(int)m,(int)col);
37647c6ae99SBarry Smith 
377aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
378aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
37947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
38047c6ae99SBarry Smith 
38147c6ae99SBarry Smith   /* create the coloring */
38247c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
38347c6ae99SBarry Smith     if (!dd->localcoloring) {
384785e854fSJed Brown       ierr = PetscMalloc1(nc*nx,&colors);CHKERRQ(ierr);
385ae4f298aSBarry Smith       if (dd->ofillcols) {
386ae4f298aSBarry Smith         PetscInt tc = 0;
387ae4f298aSBarry Smith         for (i=0; i<nc; i++) tc += (PetscInt) (dd->ofillcols[i] > 0);
388ae4f298aSBarry Smith         i1 = 0;
389ae4f298aSBarry Smith         for (i=xs; i<xs+nx; i++) {
390ae4f298aSBarry Smith           for (l=0; l<nc; l++) {
391ae4f298aSBarry Smith             if (dd->ofillcols[l] && (i % col)) {
392ae4f298aSBarry Smith               colors[i1++] =  nc - 1 + tc*((i % col) - 1) + dd->ofillcols[l];
393ae4f298aSBarry Smith             } else {
394ae4f298aSBarry Smith               colors[i1++] = l;
395ae4f298aSBarry Smith             }
396ae4f298aSBarry Smith           }
397ae4f298aSBarry Smith         }
398ae4f298aSBarry Smith         ncolors = nc + 2*s*tc;
399ae4f298aSBarry Smith       } else {
40047c6ae99SBarry Smith         i1 = 0;
40147c6ae99SBarry Smith         for (i=xs; i<xs+nx; i++) {
40247c6ae99SBarry Smith           for (l=0; l<nc; l++) {
40347c6ae99SBarry Smith             colors[i1++] = l + nc*(i % col);
40447c6ae99SBarry Smith           }
40547c6ae99SBarry Smith         }
40647c6ae99SBarry Smith         ncolors = nc + nc*(col-1);
407ae4f298aSBarry Smith       }
40847c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr);
40947c6ae99SBarry Smith     }
41047c6ae99SBarry Smith     *coloring = dd->localcoloring;
41147c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
41247c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
413785e854fSJed Brown       ierr = PetscMalloc1(nc*gnx,&colors);CHKERRQ(ierr);
41447c6ae99SBarry Smith       i1   = 0;
41547c6ae99SBarry Smith       for (i=gxs; i<gxs+gnx; i++) {
41647c6ae99SBarry Smith         for (l=0; l<nc; l++) {
41747c6ae99SBarry Smith           /* the complicated stuff is to handle periodic boundaries */
41847c6ae99SBarry Smith           colors[i1++] = l + nc*(SetInRange(i,m) % col);
41947c6ae99SBarry Smith         }
42047c6ae99SBarry Smith       }
42147c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
42247c6ae99SBarry Smith       ierr    = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
42347c6ae99SBarry Smith       ierr    = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
42447c6ae99SBarry Smith     }
42547c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
426ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
42747c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
42847c6ae99SBarry Smith   PetscFunctionReturn(0);
42947c6ae99SBarry Smith }
43047c6ae99SBarry Smith 
43147c6ae99SBarry Smith #undef __FUNCT__
432e727c939SJed Brown #define __FUNCT__ "DMCreateColoring_DA_2d_5pt_MPIAIJ"
433e727c939SJed Brown PetscErrorCode DMCreateColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
43447c6ae99SBarry Smith {
43547c6ae99SBarry Smith   PetscErrorCode   ierr;
43647c6ae99SBarry Smith   PetscInt         xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
43747c6ae99SBarry Smith   PetscInt         ncolors;
43847c6ae99SBarry Smith   MPI_Comm         comm;
4391321219cSEthan Coon   DMDABoundaryType bx,by;
44047c6ae99SBarry Smith   ISColoringValue  *colors;
44147c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
44247c6ae99SBarry Smith 
44347c6ae99SBarry Smith   PetscFunctionBegin;
44447c6ae99SBarry Smith   /*
44547c6ae99SBarry Smith          nc - number of components per grid point
44647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
44747c6ae99SBarry Smith 
44847c6ae99SBarry Smith   */
4491321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);CHKERRQ(ierr);
450aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
451aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
45247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
45347c6ae99SBarry Smith 
454ce94432eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % 5)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible by 5\n");
455ce94432eSBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % 5)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible by 5\n");
45647c6ae99SBarry Smith 
45747c6ae99SBarry Smith   /* create the coloring */
45847c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
45947c6ae99SBarry Smith     if (!dd->localcoloring) {
460785e854fSJed Brown       ierr = PetscMalloc1(nc*nx*ny,&colors);CHKERRQ(ierr);
46147c6ae99SBarry Smith       ii   = 0;
46247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
46347c6ae99SBarry Smith         for (i=xs; i<xs+nx; i++) {
46447c6ae99SBarry Smith           for (k=0; k<nc; k++) {
46547c6ae99SBarry Smith             colors[ii++] = k + nc*((3*j+i) % 5);
46647c6ae99SBarry Smith           }
46747c6ae99SBarry Smith         }
46847c6ae99SBarry Smith       }
46947c6ae99SBarry Smith       ncolors = 5*nc;
47047c6ae99SBarry Smith       ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
47147c6ae99SBarry Smith     }
47247c6ae99SBarry Smith     *coloring = dd->localcoloring;
47347c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
47447c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
475785e854fSJed Brown       ierr = PetscMalloc1(nc*gnx*gny,&colors);CHKERRQ(ierr);
47647c6ae99SBarry Smith       ii = 0;
47747c6ae99SBarry Smith       for (j=gys; j<gys+gny; j++) {
47847c6ae99SBarry Smith         for (i=gxs; i<gxs+gnx; i++) {
47947c6ae99SBarry Smith           for (k=0; k<nc; k++) {
48047c6ae99SBarry Smith             colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
48147c6ae99SBarry Smith           }
48247c6ae99SBarry Smith         }
48347c6ae99SBarry Smith       }
48447c6ae99SBarry Smith       ncolors = 5*nc;
48547c6ae99SBarry Smith       ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
48647c6ae99SBarry Smith       ierr    = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
48747c6ae99SBarry Smith     }
48847c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
489ce94432eSBarry Smith   } else SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
49047c6ae99SBarry Smith   PetscFunctionReturn(0);
49147c6ae99SBarry Smith }
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith /* =========================================================================== */
494950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM,Mat);
495ce308e1dSBarry Smith extern PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM,Mat);
496950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM,Mat);
497950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM,Mat);
498950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM,Mat);
499950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM,Mat);
500950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM,Mat);
501950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM,Mat);
502950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM,Mat);
503950540a4SJed Brown extern PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM,Mat);
50447c6ae99SBarry Smith 
50547c6ae99SBarry Smith #undef __FUNCT__
506c688c046SMatthew G Knepley #define __FUNCT__ "MatSetupDM"
5078bbdbebaSMatthew G Knepley /*@C
508c688c046SMatthew G Knepley    MatSetupDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix
50947c6ae99SBarry Smith 
51047c6ae99SBarry Smith    Logically Collective on Mat
51147c6ae99SBarry Smith 
51247c6ae99SBarry Smith    Input Parameters:
51347c6ae99SBarry Smith +  mat - the matrix
51447c6ae99SBarry Smith -  da - the da
51547c6ae99SBarry Smith 
51647c6ae99SBarry Smith    Level: intermediate
51747c6ae99SBarry Smith 
51847c6ae99SBarry Smith @*/
519c688c046SMatthew G Knepley PetscErrorCode MatSetupDM(Mat mat,DM da)
52047c6ae99SBarry Smith {
52147c6ae99SBarry Smith   PetscErrorCode ierr;
52247c6ae99SBarry Smith 
52347c6ae99SBarry Smith   PetscFunctionBegin;
52447c6ae99SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
52547c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
526c688c046SMatthew G Knepley   ierr = PetscTryMethod(mat,"MatSetupDM_C",(Mat,DM),(mat,da));CHKERRQ(ierr);
52747c6ae99SBarry Smith   PetscFunctionReturn(0);
52847c6ae99SBarry Smith }
52947c6ae99SBarry Smith 
53047c6ae99SBarry Smith #undef __FUNCT__
53147c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA"
5327087cfbeSBarry Smith PetscErrorCode  MatView_MPI_DA(Mat A,PetscViewer viewer)
53347c6ae99SBarry Smith {
5349a42bb27SBarry Smith   DM                da;
53547c6ae99SBarry Smith   PetscErrorCode    ierr;
53647c6ae99SBarry Smith   const char        *prefix;
53747c6ae99SBarry Smith   Mat               Anatural;
53847c6ae99SBarry Smith   AO                ao;
53947c6ae99SBarry Smith   PetscInt          rstart,rend,*petsc,i;
54047c6ae99SBarry Smith   IS                is;
54147c6ae99SBarry Smith   MPI_Comm          comm;
54274388724SJed Brown   PetscViewerFormat format;
54347c6ae99SBarry Smith 
54447c6ae99SBarry Smith   PetscFunctionBegin;
54574388724SJed Brown   /* Check whether we are just printing info, in which case MatView() already viewed everything we wanted to view */
54674388724SJed Brown   ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr);
54774388724SJed Brown   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) PetscFunctionReturn(0);
54874388724SJed Brown 
54947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
550c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
551ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
55247c6ae99SBarry Smith 
553aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
55447c6ae99SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
555785e854fSJed Brown   ierr = PetscMalloc1((rend-rstart),&petsc);CHKERRQ(ierr);
55647c6ae99SBarry Smith   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
55747c6ae99SBarry Smith   ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr);
55847c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
55947c6ae99SBarry Smith 
56047c6ae99SBarry Smith   /* call viewer on natural ordering */
56147c6ae99SBarry Smith   ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr);
562fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
56347c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr);
56447c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr);
56547c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr);
56647c6ae99SBarry Smith   ierr = MatView(Anatural,viewer);CHKERRQ(ierr);
567fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
56847c6ae99SBarry Smith   PetscFunctionReturn(0);
56947c6ae99SBarry Smith }
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith #undef __FUNCT__
57247c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA"
5737087cfbeSBarry Smith PetscErrorCode  MatLoad_MPI_DA(Mat A,PetscViewer viewer)
57447c6ae99SBarry Smith {
5759a42bb27SBarry Smith   DM             da;
57647c6ae99SBarry Smith   PetscErrorCode ierr;
57747c6ae99SBarry Smith   Mat            Anatural,Aapp;
57847c6ae99SBarry Smith   AO             ao;
57947c6ae99SBarry Smith   PetscInt       rstart,rend,*app,i;
58047c6ae99SBarry Smith   IS             is;
58147c6ae99SBarry Smith   MPI_Comm       comm;
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith   PetscFunctionBegin;
58447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
585c688c046SMatthew G Knepley   ierr = MatGetDM(A, &da);CHKERRQ(ierr);
586ce94432eSBarry Smith   if (!da) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
58747c6ae99SBarry Smith 
58847c6ae99SBarry Smith   /* Load the matrix in natural ordering */
589ce94432eSBarry Smith   ierr = MatCreate(PetscObjectComm((PetscObject)A),&Anatural);CHKERRQ(ierr);
59047c6ae99SBarry Smith   ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr);
59147c6ae99SBarry Smith   ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
59247c6ae99SBarry Smith   ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr);
59347c6ae99SBarry Smith 
59447c6ae99SBarry Smith   /* Map natural ordering to application ordering and create IS */
595aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
59647c6ae99SBarry Smith   ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr);
597785e854fSJed Brown   ierr = PetscMalloc1((rend-rstart),&app);CHKERRQ(ierr);
59847c6ae99SBarry Smith   for (i=rstart; i<rend; i++) app[i-rstart] = i;
59947c6ae99SBarry Smith   ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr);
60047c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
60147c6ae99SBarry Smith 
60247c6ae99SBarry Smith   /* Do permutation and replace header */
60347c6ae99SBarry Smith   ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr);
60447c6ae99SBarry Smith   ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr);
605fcfd50ebSBarry Smith   ierr = ISDestroy(&is);CHKERRQ(ierr);
606fcfd50ebSBarry Smith   ierr = MatDestroy(&Anatural);CHKERRQ(ierr);
60747c6ae99SBarry Smith   PetscFunctionReturn(0);
60847c6ae99SBarry Smith }
60947c6ae99SBarry Smith 
61047c6ae99SBarry Smith #undef __FUNCT__
611950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA"
612b412c318SBarry Smith PetscErrorCode DMCreateMatrix_DA(DM da, Mat *J)
61347c6ae99SBarry Smith {
61447c6ae99SBarry Smith   PetscErrorCode ierr;
61547c6ae99SBarry Smith   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
61647c6ae99SBarry Smith   Mat            A;
61747c6ae99SBarry Smith   MPI_Comm       comm;
61819fd82e9SBarry Smith   MatType        Atype;
61937d0c07bSMatthew G Knepley   PetscSection   section, sectionGlobal;
6200298fd71SBarry Smith   void           (*aij)(void)=NULL,(*baij)(void)=NULL,(*sbaij)(void)=NULL;
621b412c318SBarry Smith   MatType        mtype;
62247c6ae99SBarry Smith   PetscMPIInt    size;
62347c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
62447c6ae99SBarry Smith 
62547c6ae99SBarry Smith   PetscFunctionBegin;
626607a6623SBarry Smith   ierr = MatInitializePackage();CHKERRQ(ierr);
627b412c318SBarry Smith   mtype = da->mattype;
62847c6ae99SBarry Smith 
62937d0c07bSMatthew G Knepley   ierr = DMGetDefaultSection(da, &section);CHKERRQ(ierr);
63037d0c07bSMatthew G Knepley   if (section) {
63137d0c07bSMatthew G Knepley     PetscInt  bs = -1;
63237d0c07bSMatthew G Knepley     PetscInt  localSize;
63337d0c07bSMatthew G Knepley     PetscBool isShell, isBlock, isSeqBlock, isMPIBlock, isSymBlock, isSymSeqBlock, isSymMPIBlock, isSymmetric;
63437d0c07bSMatthew G Knepley 
63537d0c07bSMatthew G Knepley     ierr = DMGetDefaultGlobalSection(da, &sectionGlobal);CHKERRQ(ierr);
63637d0c07bSMatthew G Knepley     ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal, &localSize);CHKERRQ(ierr);
63782f516ccSBarry Smith     ierr = MatCreate(PetscObjectComm((PetscObject)da), J);CHKERRQ(ierr);
63837d0c07bSMatthew G Knepley     ierr = MatSetSizes(*J, localSize, localSize, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
63937d0c07bSMatthew G Knepley     ierr = MatSetType(*J, mtype);CHKERRQ(ierr);
64037d0c07bSMatthew G Knepley     ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
64137d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSHELL, &isShell);CHKERRQ(ierr);
64237d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATBAIJ, &isBlock);CHKERRQ(ierr);
64337d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSEQBAIJ, &isSeqBlock);CHKERRQ(ierr);
64437d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATMPIBAIJ, &isMPIBlock);CHKERRQ(ierr);
64537d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr);
64637d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr);
64737d0c07bSMatthew G Knepley     ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr);
64837d0c07bSMatthew G Knepley     /* Check for symmetric storage */
64937d0c07bSMatthew G Knepley     isSymmetric = (PetscBool) (isSymBlock || isSymSeqBlock || isSymMPIBlock);
65037d0c07bSMatthew G Knepley     if (isSymmetric) {
65137d0c07bSMatthew G Knepley       ierr = MatSetOption(*J, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);
65237d0c07bSMatthew G Knepley     }
65337d0c07bSMatthew G Knepley     if (!isShell) {
65437d0c07bSMatthew G Knepley       PetscInt *dnz, *onz, *dnzu, *onzu, bsLocal;
65537d0c07bSMatthew G Knepley 
65637d0c07bSMatthew G Knepley       if (bs < 0) {
65737d0c07bSMatthew G Knepley         if (isBlock || isSeqBlock || isMPIBlock || isSymBlock || isSymSeqBlock || isSymMPIBlock) {
65837d0c07bSMatthew G Knepley           PetscInt pStart, pEnd, p, dof;
65937d0c07bSMatthew G Knepley 
66037d0c07bSMatthew G Knepley           ierr = PetscSectionGetChart(sectionGlobal, &pStart, &pEnd);CHKERRQ(ierr);
66137d0c07bSMatthew G Knepley           for (p = pStart; p < pEnd; ++p) {
66237d0c07bSMatthew G Knepley             ierr = PetscSectionGetDof(sectionGlobal, p, &dof);CHKERRQ(ierr);
66337d0c07bSMatthew G Knepley             if (dof) {
66437d0c07bSMatthew G Knepley               bs = dof;
66537d0c07bSMatthew G Knepley               break;
66637d0c07bSMatthew G Knepley             }
66737d0c07bSMatthew G Knepley           }
66837d0c07bSMatthew G Knepley         } else {
66937d0c07bSMatthew G Knepley           bs = 1;
67037d0c07bSMatthew G Knepley         }
67137d0c07bSMatthew G Knepley         /* Must have same blocksize on all procs (some might have no points) */
67237d0c07bSMatthew G Knepley         bsLocal = bs;
67382f516ccSBarry Smith         ierr    = MPI_Allreduce(&bsLocal, &bs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
67437d0c07bSMatthew G Knepley       }
675*1795a4d1SJed Brown       ierr = PetscCalloc4(localSize/bs, &dnz, localSize/bs, &onz, localSize/bs, &dnzu, localSize/bs, &onzu);CHKERRQ(ierr);
676552f7358SJed Brown       /* ierr = DMPlexPreallocateOperator(dm, bs, section, sectionGlobal, dnz, onz, dnzu, onzu, *J, fillMatrix);CHKERRQ(ierr); */
67737d0c07bSMatthew G Knepley       ierr = PetscFree4(dnz, onz, dnzu, onzu);CHKERRQ(ierr);
67837d0c07bSMatthew G Knepley     }
67937d0c07bSMatthew G Knepley   }
68047c6ae99SBarry Smith   /*
68147c6ae99SBarry Smith                                   m
68247c6ae99SBarry Smith           ------------------------------------------------------
68347c6ae99SBarry Smith          |                                                     |
68447c6ae99SBarry Smith          |                                                     |
68547c6ae99SBarry Smith          |               ----------------------                |
68647c6ae99SBarry Smith          |               |                    |                |
68747c6ae99SBarry Smith       n  |           ny  |                    |                |
68847c6ae99SBarry Smith          |               |                    |                |
68947c6ae99SBarry Smith          |               .---------------------                |
69047c6ae99SBarry Smith          |             (xs,ys)     nx                          |
69147c6ae99SBarry Smith          |            .                                        |
69247c6ae99SBarry Smith          |         (gxs,gys)                                   |
69347c6ae99SBarry Smith          |                                                     |
69447c6ae99SBarry Smith           -----------------------------------------------------
69547c6ae99SBarry Smith   */
69647c6ae99SBarry Smith 
69747c6ae99SBarry Smith   /*
69847c6ae99SBarry Smith          nc - number of components per grid point
69947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
70047c6ae99SBarry Smith 
70147c6ae99SBarry Smith   */
702e30e807fSPeter Brune   M   = dd->M;
703e30e807fSPeter Brune   N   = dd->N;
704e30e807fSPeter Brune   P   = dd->P;
705e30e807fSPeter Brune   dim = dd->dim;
706e30e807fSPeter Brune   dof = dd->w;
707e30e807fSPeter Brune   /* ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); */
708aa219208SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
70947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
71047c6ae99SBarry Smith   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
71147c6ae99SBarry Smith   ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr);
712b412c318SBarry Smith   ierr = MatSetType(A,mtype);CHKERRQ(ierr);
71395ee5b0eSBarry Smith   ierr = MatSetDM(A,da);CHKERRQ(ierr);
71447c6ae99SBarry Smith   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
71547c6ae99SBarry Smith   ierr = MatGetType(A,&Atype);CHKERRQ(ierr);
71647c6ae99SBarry Smith   /*
717aa219208SBarry Smith      We do not provide a getmatrix function in the DMDA operations because
718aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
71947c6ae99SBarry Smith    more low-level than matrices. This is kind of cheating but, cause sometimes
720aa219208SBarry Smith    we think of DMDA has higher level than matrices.
72147c6ae99SBarry Smith 
72247c6ae99SBarry Smith      We could switch based on Atype (or mtype), but we do not since the
72347c6ae99SBarry Smith    specialized setting routines depend only the particular preallocation
72447c6ae99SBarry Smith    details of the matrix, not the type itself.
72547c6ae99SBarry Smith   */
72647c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
72747c6ae99SBarry Smith   if (!aij) {
72847c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
72947c6ae99SBarry Smith   }
73047c6ae99SBarry Smith   if (!aij) {
73147c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
73247c6ae99SBarry Smith     if (!baij) {
73347c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
73447c6ae99SBarry Smith     }
73547c6ae99SBarry Smith     if (!baij) {
73647c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
73747c6ae99SBarry Smith       if (!sbaij) {
73847c6ae99SBarry Smith         ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
73947c6ae99SBarry Smith       }
74047c6ae99SBarry Smith     }
74147c6ae99SBarry Smith   }
74247c6ae99SBarry Smith   if (aij) {
74347c6ae99SBarry Smith     if (dim == 1) {
744ce308e1dSBarry Smith       if (dd->ofill) {
745ce308e1dSBarry Smith         ierr = DMCreateMatrix_DA_1d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
746ce308e1dSBarry Smith       } else {
747950540a4SJed Brown         ierr = DMCreateMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr);
748ce308e1dSBarry Smith       }
74947c6ae99SBarry Smith     } else if (dim == 2) {
75047c6ae99SBarry Smith       if (dd->ofill) {
751950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
75247c6ae99SBarry Smith       } else {
753950540a4SJed Brown         ierr = DMCreateMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr);
75447c6ae99SBarry Smith       }
75547c6ae99SBarry Smith     } else if (dim == 3) {
75647c6ae99SBarry Smith       if (dd->ofill) {
757950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
75847c6ae99SBarry Smith       } else {
759950540a4SJed Brown         ierr = DMCreateMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr);
76047c6ae99SBarry Smith       }
76147c6ae99SBarry Smith     }
76247c6ae99SBarry Smith   } else if (baij) {
76347c6ae99SBarry Smith     if (dim == 2) {
764950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr);
76547c6ae99SBarry Smith     } else if (dim == 3) {
766950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr);
767ce94432eSBarry 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);
76847c6ae99SBarry Smith   } else if (sbaij) {
76947c6ae99SBarry Smith     if (dim == 2) {
770950540a4SJed Brown       ierr = DMCreateMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr);
77147c6ae99SBarry Smith     } else if (dim == 3) {
772950540a4SJed Brown       ierr = DMCreateMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr);
773ce94432eSBarry 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);
774869776cdSLisandro Dalcin   } else {
775869776cdSLisandro Dalcin     ISLocalToGlobalMapping ltog,ltogb;
776869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
777869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
7782949035bSJed Brown     ierr = MatSetUp(A);CHKERRQ(ierr);
779869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMapping(A,ltog,ltog);CHKERRQ(ierr);
780869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMappingBlock(A,ltogb,ltogb);CHKERRQ(ierr);
78147c6ae99SBarry Smith   }
782aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
78347c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
784c688c046SMatthew G Knepley   ierr = MatSetDM(A,da);CHKERRQ(ierr);
78547c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
78647c6ae99SBarry Smith   if (size > 1) {
78747c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
78847c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MPI_DA);CHKERRQ(ierr);
78947c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void))MatLoad_MPI_DA);CHKERRQ(ierr);
79047c6ae99SBarry Smith   }
79147c6ae99SBarry Smith   *J = A;
79247c6ae99SBarry Smith   PetscFunctionReturn(0);
79347c6ae99SBarry Smith }
79447c6ae99SBarry Smith 
79547c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
79647c6ae99SBarry Smith #undef __FUNCT__
797950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ"
798950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM da,Mat J)
79947c6ae99SBarry Smith {
80047c6ae99SBarry Smith   PetscErrorCode         ierr;
8010298fd71SBarry 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;
80247c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
80347c6ae99SBarry Smith   MPI_Comm               comm;
80447c6ae99SBarry Smith   PetscScalar            *values;
8051321219cSEthan Coon   DMDABoundaryType       bx,by;
80647c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
807aa219208SBarry Smith   DMDAStencilType        st;
80847c6ae99SBarry Smith 
80947c6ae99SBarry Smith   PetscFunctionBegin;
81047c6ae99SBarry Smith   /*
81147c6ae99SBarry Smith          nc - number of components per grid point
81247c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
81347c6ae99SBarry Smith 
81447c6ae99SBarry Smith   */
8151321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
81647c6ae99SBarry Smith   col  = 2*s + 1;
817aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
818aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
81947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
82047c6ae99SBarry Smith 
821dcca6d9dSJed Brown   ierr = PetscMalloc2(nc,&rows,col*col*nc*nc,&cols);CHKERRQ(ierr);
8221411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
8231411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
82447c6ae99SBarry Smith 
82547c6ae99SBarry Smith   /* determine the matrix preallocation information */
82647c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
82747c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
82847c6ae99SBarry Smith 
8291321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8301321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
83147c6ae99SBarry Smith 
83247c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
83347c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
83447c6ae99SBarry Smith 
8351321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8361321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
83747c6ae99SBarry Smith 
83847c6ae99SBarry Smith       cnt = 0;
83947c6ae99SBarry Smith       for (k=0; k<nc; k++) {
84047c6ae99SBarry Smith         for (l=lstart; l<lend+1; l++) {
84147c6ae99SBarry Smith           for (p=pstart; p<pend+1; p++) {
842aa219208SBarry Smith             if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
84347c6ae99SBarry Smith               cols[cnt++] = k + nc*(slot + gnx*l + p);
84447c6ae99SBarry Smith             }
84547c6ae99SBarry Smith           }
84647c6ae99SBarry Smith         }
84747c6ae99SBarry Smith         rows[k] = k + nc*(slot);
84847c6ae99SBarry Smith       }
849784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
85047c6ae99SBarry Smith     }
85147c6ae99SBarry Smith   }
852f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
85347c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
85447c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
85547c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
85647c6ae99SBarry Smith 
857784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
858784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
85947c6ae99SBarry Smith 
86047c6ae99SBarry Smith   /*
86147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
86247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
86347c6ae99SBarry Smith     PETSc ordering.
86447c6ae99SBarry Smith   */
865fcfd50ebSBarry Smith   if (!da->prealloc_only) {
866*1795a4d1SJed Brown     ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr);
86747c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
86847c6ae99SBarry Smith 
8691321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8701321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
87147c6ae99SBarry Smith 
87247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
87347c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
87447c6ae99SBarry Smith 
8751321219cSEthan Coon         lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8761321219cSEthan Coon         lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
87747c6ae99SBarry Smith 
87847c6ae99SBarry Smith         cnt = 0;
87947c6ae99SBarry Smith         for (k=0; k<nc; k++) {
88047c6ae99SBarry Smith           for (l=lstart; l<lend+1; l++) {
88147c6ae99SBarry Smith             for (p=pstart; p<pend+1; p++) {
882aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
88347c6ae99SBarry Smith                 cols[cnt++] = k + nc*(slot + gnx*l + p);
88447c6ae99SBarry Smith               }
88547c6ae99SBarry Smith             }
88647c6ae99SBarry Smith           }
88747c6ae99SBarry Smith           rows[k] = k + nc*(slot);
88847c6ae99SBarry Smith         }
88947c6ae99SBarry Smith         ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
89047c6ae99SBarry Smith       }
89147c6ae99SBarry Smith     }
89247c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
89347c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89447c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89547c6ae99SBarry Smith   }
89647c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
89747c6ae99SBarry Smith   PetscFunctionReturn(0);
89847c6ae99SBarry Smith }
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith #undef __FUNCT__
901950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ_Fill"
902950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J)
90347c6ae99SBarry Smith {
90447c6ae99SBarry Smith   PetscErrorCode         ierr;
90547c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
90647c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p;
90747c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
90847c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
90947c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
91047c6ae99SBarry Smith   MPI_Comm               comm;
91147c6ae99SBarry Smith   PetscScalar            *values;
9121321219cSEthan Coon   DMDABoundaryType       bx,by;
91347c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
914aa219208SBarry Smith   DMDAStencilType        st;
91547c6ae99SBarry Smith 
91647c6ae99SBarry Smith   PetscFunctionBegin;
91747c6ae99SBarry Smith   /*
91847c6ae99SBarry Smith          nc - number of components per grid point
91947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
92047c6ae99SBarry Smith 
92147c6ae99SBarry Smith   */
9221321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
92347c6ae99SBarry Smith   col  = 2*s + 1;
924aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
925aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
92647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
92747c6ae99SBarry Smith 
928785e854fSJed Brown   ierr = PetscMalloc1(col*col*nc*nc,&cols);CHKERRQ(ierr);
9291411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
9301411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
93147c6ae99SBarry Smith 
93247c6ae99SBarry Smith   /* determine the matrix preallocation information */
93347c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
93447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
93547c6ae99SBarry Smith 
9361321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9371321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
93847c6ae99SBarry Smith 
93947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
94047c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
94147c6ae99SBarry Smith 
9421321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9431321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
94447c6ae99SBarry Smith 
94547c6ae99SBarry Smith       for (k=0; k<nc; k++) {
94647c6ae99SBarry Smith         cnt = 0;
94747c6ae99SBarry Smith         for (l=lstart; l<lend+1; l++) {
94847c6ae99SBarry Smith           for (p=pstart; p<pend+1; p++) {
94947c6ae99SBarry Smith             if (l || p) {
950aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
9518865f1eaSKarl Rupp                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p);
95247c6ae99SBarry Smith               }
95347c6ae99SBarry Smith             } else {
95447c6ae99SBarry Smith               if (dfill) {
9558865f1eaSKarl Rupp                 for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p);
95647c6ae99SBarry Smith               } else {
9578865f1eaSKarl Rupp                 for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + gnx*l + p);
95847c6ae99SBarry Smith               }
95947c6ae99SBarry Smith             }
96047c6ae99SBarry Smith           }
96147c6ae99SBarry Smith         }
96247c6ae99SBarry Smith         row  = k + nc*(slot);
963784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
96447c6ae99SBarry Smith       }
96547c6ae99SBarry Smith     }
96647c6ae99SBarry Smith   }
96747c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
96847c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
96947c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
970784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
971784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
97247c6ae99SBarry Smith 
97347c6ae99SBarry Smith   /*
97447c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
97547c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
97647c6ae99SBarry Smith     PETSc ordering.
97747c6ae99SBarry Smith   */
978fcfd50ebSBarry Smith   if (!da->prealloc_only) {
979*1795a4d1SJed Brown     ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr);
98047c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
98147c6ae99SBarry Smith 
9821321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9831321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
98447c6ae99SBarry Smith 
98547c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
98647c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
98747c6ae99SBarry Smith 
9881321219cSEthan Coon         lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9891321219cSEthan Coon         lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
99047c6ae99SBarry Smith 
99147c6ae99SBarry Smith         for (k=0; k<nc; k++) {
99247c6ae99SBarry Smith           cnt = 0;
99347c6ae99SBarry Smith           for (l=lstart; l<lend+1; l++) {
99447c6ae99SBarry Smith             for (p=pstart; p<pend+1; p++) {
99547c6ae99SBarry Smith               if (l || p) {
996aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
9978865f1eaSKarl Rupp                   for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p);
99847c6ae99SBarry Smith                 }
99947c6ae99SBarry Smith               } else {
100047c6ae99SBarry Smith                 if (dfill) {
10018865f1eaSKarl Rupp                   for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p);
100247c6ae99SBarry Smith                 } else {
10038865f1eaSKarl Rupp                   for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + gnx*l + p);
100447c6ae99SBarry Smith                 }
100547c6ae99SBarry Smith               }
100647c6ae99SBarry Smith             }
100747c6ae99SBarry Smith           }
100847c6ae99SBarry Smith           row  = k + nc*(slot);
100947c6ae99SBarry Smith           ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
101047c6ae99SBarry Smith         }
101147c6ae99SBarry Smith       }
101247c6ae99SBarry Smith     }
101347c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
101447c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101547c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101647c6ae99SBarry Smith   }
101747c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
101847c6ae99SBarry Smith   PetscFunctionReturn(0);
101947c6ae99SBarry Smith }
102047c6ae99SBarry Smith 
102147c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
102247c6ae99SBarry Smith 
102347c6ae99SBarry Smith #undef __FUNCT__
1024950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ"
1025950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM da,Mat J)
102647c6ae99SBarry Smith {
102747c6ae99SBarry Smith   PetscErrorCode         ierr;
102847c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
10290298fd71SBarry Smith   PetscInt               m,n,dim,s,*cols = NULL,k,nc,*rows = NULL,col,cnt,l,p,*dnz = NULL,*onz = NULL;
103047c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
103147c6ae99SBarry Smith   MPI_Comm               comm;
103247c6ae99SBarry Smith   PetscScalar            *values;
10331321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
103447c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1035aa219208SBarry Smith   DMDAStencilType        st;
103647c6ae99SBarry Smith 
103747c6ae99SBarry Smith   PetscFunctionBegin;
103847c6ae99SBarry Smith   /*
103947c6ae99SBarry Smith          nc - number of components per grid point
104047c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
104147c6ae99SBarry Smith 
104247c6ae99SBarry Smith   */
10431321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
104447c6ae99SBarry Smith   col  = 2*s + 1;
104547c6ae99SBarry Smith 
1046aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1047aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
104847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
104947c6ae99SBarry Smith 
1050dcca6d9dSJed Brown   ierr = PetscMalloc2(nc,&rows,col*col*col*nc*nc,&cols);CHKERRQ(ierr);
10511411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
10521411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
105347c6ae99SBarry Smith 
105447c6ae99SBarry Smith   /* determine the matrix preallocation information */
105547c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
105647c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
10571321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10581321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
105947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
10601321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10611321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
106247c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
10631321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
10641321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
106547c6ae99SBarry Smith 
106647c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
106747c6ae99SBarry Smith 
106847c6ae99SBarry Smith         cnt = 0;
106947c6ae99SBarry Smith         for (l=0; l<nc; l++) {
107047c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
107147c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
107247c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1073aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
107447c6ae99SBarry Smith                   cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
107547c6ae99SBarry Smith                 }
107647c6ae99SBarry Smith               }
107747c6ae99SBarry Smith             }
107847c6ae99SBarry Smith           }
107947c6ae99SBarry Smith           rows[l] = l + nc*(slot);
108047c6ae99SBarry Smith         }
1081784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
108247c6ae99SBarry Smith       }
108347c6ae99SBarry Smith     }
108447c6ae99SBarry Smith   }
1085f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
108647c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
108747c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
108847c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1089784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1090784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
109147c6ae99SBarry Smith 
109247c6ae99SBarry Smith   /*
109347c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
109447c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
109547c6ae99SBarry Smith     PETSc ordering.
109647c6ae99SBarry Smith   */
1097fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1098*1795a4d1SJed Brown     ierr = PetscCalloc1(col*col*col*nc*nc*nc,&values);CHKERRQ(ierr);
109947c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
11001321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
11011321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
110247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
11031321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
11041321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
110547c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
11061321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
11071321219cSEthan Coon           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
110847c6ae99SBarry Smith 
110947c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
111047c6ae99SBarry Smith 
111147c6ae99SBarry Smith           cnt = 0;
111247c6ae99SBarry Smith           for (l=0; l<nc; l++) {
111347c6ae99SBarry Smith             for (ii=istart; ii<iend+1; ii++) {
111447c6ae99SBarry Smith               for (jj=jstart; jj<jend+1; jj++) {
111547c6ae99SBarry Smith                 for (kk=kstart; kk<kend+1; kk++) {
1116aa219208SBarry Smith                   if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
111747c6ae99SBarry Smith                     cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
111847c6ae99SBarry Smith                   }
111947c6ae99SBarry Smith                 }
112047c6ae99SBarry Smith               }
112147c6ae99SBarry Smith             }
112247c6ae99SBarry Smith             rows[l] = l + nc*(slot);
112347c6ae99SBarry Smith           }
112447c6ae99SBarry Smith           ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
112547c6ae99SBarry Smith         }
112647c6ae99SBarry Smith       }
112747c6ae99SBarry Smith     }
112847c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
112947c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
113047c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
113147c6ae99SBarry Smith   }
113247c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
113347c6ae99SBarry Smith   PetscFunctionReturn(0);
113447c6ae99SBarry Smith }
113547c6ae99SBarry Smith 
113647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
113747c6ae99SBarry Smith 
113847c6ae99SBarry Smith #undef __FUNCT__
1139ce308e1dSBarry Smith #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ_Fill"
1140ce308e1dSBarry Smith PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM da,Mat J)
1141ce308e1dSBarry Smith {
1142ce308e1dSBarry Smith   PetscErrorCode         ierr;
1143ce308e1dSBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1144ce308e1dSBarry Smith   PetscInt               xs,nx,i,j,gxs,gnx,row,k,l;
11450298fd71SBarry Smith   PetscInt               m,dim,s,*cols = NULL,nc,col,cnt,*ocols;
1146ce308e1dSBarry Smith   PetscInt               *ofill = dd->ofill;
1147ce308e1dSBarry Smith   PetscScalar            *values;
1148ce308e1dSBarry Smith   DMDABoundaryType       bx;
1149ce308e1dSBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1150ce308e1dSBarry Smith   PetscMPIInt            rank,size;
1151ce308e1dSBarry Smith 
1152ce308e1dSBarry Smith   PetscFunctionBegin;
1153ce94432eSBarry Smith   if (dd->bx == DMDA_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"With fill provided not implemented with periodic boundary conditions");
1154ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
1155ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr);
1156ce308e1dSBarry Smith 
1157ce308e1dSBarry Smith   /*
1158ce308e1dSBarry Smith          nc - number of components per grid point
1159ce308e1dSBarry Smith          col - number of colors needed in one direction for single component problem
1160ce308e1dSBarry Smith 
1161ce308e1dSBarry Smith   */
1162ce308e1dSBarry Smith   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
1163ce308e1dSBarry Smith   col  = 2*s + 1;
1164ce308e1dSBarry Smith 
1165ce308e1dSBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1166ce308e1dSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
1167ce308e1dSBarry Smith 
1168ce308e1dSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
1169*1795a4d1SJed Brown   ierr = PetscCalloc2(nx*nc,&cols,nx*nc,&ocols);CHKERRQ(ierr);
1170ce308e1dSBarry Smith 
1171ce308e1dSBarry Smith   /*
1172ce308e1dSBarry Smith         note should be smaller for first and last process with no periodic
1173ce308e1dSBarry Smith         does not handle dfill
1174ce308e1dSBarry Smith   */
1175ce308e1dSBarry Smith   cnt = 0;
1176ce308e1dSBarry Smith   /* coupling with process to the left */
1177ce308e1dSBarry Smith   for (i=0; i<s; i++) {
1178ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
1179ce308e1dSBarry Smith       ocols[cnt] = ((!rank) ? 0 : (s - i)*(ofill[j+1] - ofill[j]));
1180ce308e1dSBarry Smith       cols[cnt]  = nc + (s + i)*(ofill[j+1] - ofill[j]);
1181ce308e1dSBarry Smith       cnt++;
1182ce308e1dSBarry Smith     }
1183ce308e1dSBarry Smith   }
1184ce308e1dSBarry Smith   for (i=s; i<nx-s; i++) {
1185ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
1186ce308e1dSBarry Smith       cols[cnt] = nc + 2*s*(ofill[j+1] - ofill[j]);
1187ce308e1dSBarry Smith       cnt++;
1188ce308e1dSBarry Smith     }
1189ce308e1dSBarry Smith   }
1190ce308e1dSBarry Smith   /* coupling with process to the right */
1191ce308e1dSBarry Smith   for (i=nx-s; i<nx; i++) {
1192ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
1193ce308e1dSBarry Smith       ocols[cnt] = ((rank == (size-1)) ? 0 : (i - nx + s + 1)*(ofill[j+1] - ofill[j]));
1194ce308e1dSBarry Smith       cols[cnt]  = nc + (s + nx - i - 1)*(ofill[j+1] - ofill[j]);
1195ce308e1dSBarry Smith       cnt++;
1196ce308e1dSBarry Smith     }
1197ce308e1dSBarry Smith   }
1198ce308e1dSBarry Smith 
1199ce308e1dSBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,cols);CHKERRQ(ierr);
1200ce308e1dSBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,cols,0,ocols);CHKERRQ(ierr);
1201ce308e1dSBarry Smith   ierr = PetscFree2(cols,ocols);CHKERRQ(ierr);
1202ce308e1dSBarry Smith 
1203ce308e1dSBarry Smith   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1204ce308e1dSBarry Smith   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1205ce308e1dSBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1206ce308e1dSBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
1207ce308e1dSBarry Smith 
1208ce308e1dSBarry Smith   /*
1209ce308e1dSBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1210ce308e1dSBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1211ce308e1dSBarry Smith     PETSc ordering.
1212ce308e1dSBarry Smith   */
1213ce308e1dSBarry Smith   if (!da->prealloc_only) {
1214785e854fSJed Brown     ierr = PetscMalloc1(col*nc*nc,&cols);CHKERRQ(ierr);
1215*1795a4d1SJed Brown     ierr = PetscCalloc1(col*nc*nc,&values);CHKERRQ(ierr);
1216ce308e1dSBarry Smith 
1217ce308e1dSBarry Smith     row = xs*nc;
1218ce308e1dSBarry Smith     /* coupling with process to the left */
1219ce308e1dSBarry Smith     for (i=xs; i<xs+s; i++) {
1220ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1221ce308e1dSBarry Smith         cnt = 0;
1222ce308e1dSBarry Smith         if (rank) {
1223ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1224ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1225ce308e1dSBarry Smith           }
1226ce308e1dSBarry Smith         }
1227ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1228ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1229ce308e1dSBarry Smith         }
1230ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1231ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1232ce308e1dSBarry Smith         }
1233ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1234ce308e1dSBarry Smith         row++;
1235ce308e1dSBarry Smith       }
1236ce308e1dSBarry Smith     }
1237ce308e1dSBarry Smith     for (i=xs+s; i<xs+nx-s; i++) {
1238ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1239ce308e1dSBarry Smith         cnt = 0;
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         for (k=0; k<nc; k++) {
1244ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1245ce308e1dSBarry Smith         }
1246ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1247ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1248ce308e1dSBarry Smith         }
1249ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1250ce308e1dSBarry Smith         row++;
1251ce308e1dSBarry Smith       }
1252ce308e1dSBarry Smith     }
1253ce308e1dSBarry Smith     /* coupling with process to the right */
1254ce308e1dSBarry Smith     for (i=xs+nx-s; i<xs+nx; i++) {
1255ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1256ce308e1dSBarry Smith         cnt = 0;
1257ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1258ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1259ce308e1dSBarry Smith         }
1260ce308e1dSBarry Smith         for (k=0; k<nc; k++) {
1261ce308e1dSBarry Smith           cols[cnt++] = i*nc + k;
1262ce308e1dSBarry Smith         }
1263ce308e1dSBarry Smith         if (rank < size-1) {
1264ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1265ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1266ce308e1dSBarry Smith           }
1267ce308e1dSBarry Smith         }
1268ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1269ce308e1dSBarry Smith         row++;
1270ce308e1dSBarry Smith       }
1271ce308e1dSBarry Smith     }
1272ce308e1dSBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1273ce308e1dSBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1274ce308e1dSBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1275ce308e1dSBarry Smith     ierr = PetscFree(cols);CHKERRQ(ierr);
1276ce308e1dSBarry Smith   }
1277ce308e1dSBarry Smith   PetscFunctionReturn(0);
1278ce308e1dSBarry Smith }
1279ce308e1dSBarry Smith 
1280ce308e1dSBarry Smith /* ---------------------------------------------------------------------------------*/
1281ce308e1dSBarry Smith 
1282ce308e1dSBarry Smith #undef __FUNCT__
1283950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ"
1284950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM da,Mat J)
128547c6ae99SBarry Smith {
128647c6ae99SBarry Smith   PetscErrorCode         ierr;
128747c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
12880298fd71SBarry Smith   PetscInt               m,dim,s,*cols = NULL,nc,*rows = NULL,col,cnt,l;
128947c6ae99SBarry Smith   PetscInt               istart,iend;
129047c6ae99SBarry Smith   PetscScalar            *values;
12911321219cSEthan Coon   DMDABoundaryType       bx;
129247c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
129347c6ae99SBarry Smith 
129447c6ae99SBarry Smith   PetscFunctionBegin;
129547c6ae99SBarry Smith   /*
129647c6ae99SBarry Smith          nc - number of components per grid point
129747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
129847c6ae99SBarry Smith 
129947c6ae99SBarry Smith   */
13001321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
130147c6ae99SBarry Smith   col  = 2*s + 1;
130247c6ae99SBarry Smith 
1303aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1304aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
130547c6ae99SBarry Smith 
1306f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
130747c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
130847c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
130947c6ae99SBarry Smith 
13101411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
13111411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1312784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1313784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
131447c6ae99SBarry Smith 
131547c6ae99SBarry Smith   /*
131647c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
131747c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
131847c6ae99SBarry Smith     PETSc ordering.
131947c6ae99SBarry Smith   */
1320fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1321dcca6d9dSJed Brown     ierr = PetscMalloc2(nc,&rows,col*nc*nc,&cols);CHKERRQ(ierr);
1322*1795a4d1SJed Brown     ierr = PetscCalloc1(col*nc*nc,&values);CHKERRQ(ierr);
132347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
132447c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
132547c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
132647c6ae99SBarry Smith       slot   = i - gxs;
132747c6ae99SBarry Smith 
132847c6ae99SBarry Smith       cnt = 0;
132947c6ae99SBarry Smith       for (l=0; l<nc; l++) {
133047c6ae99SBarry Smith         for (i1=istart; i1<iend+1; i1++) {
133147c6ae99SBarry Smith           cols[cnt++] = l + nc*(slot + i1);
133247c6ae99SBarry Smith         }
133347c6ae99SBarry Smith         rows[l] = l + nc*(slot);
133447c6ae99SBarry Smith       }
133547c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
133647c6ae99SBarry Smith     }
133747c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
133847c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
133947c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
134047c6ae99SBarry Smith     ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
1341ce308e1dSBarry Smith   }
134247c6ae99SBarry Smith   PetscFunctionReturn(0);
134347c6ae99SBarry Smith }
134447c6ae99SBarry Smith 
134547c6ae99SBarry Smith #undef __FUNCT__
1346950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIBAIJ"
1347950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM da,Mat J)
134847c6ae99SBarry Smith {
134947c6ae99SBarry Smith   PetscErrorCode         ierr;
135047c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
135147c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
135247c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
135347c6ae99SBarry Smith   MPI_Comm               comm;
135447c6ae99SBarry Smith   PetscScalar            *values;
13551321219cSEthan Coon   DMDABoundaryType       bx,by;
1356aa219208SBarry Smith   DMDAStencilType        st;
135747c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
135847c6ae99SBarry Smith 
135947c6ae99SBarry Smith   PetscFunctionBegin;
136047c6ae99SBarry Smith   /*
136147c6ae99SBarry Smith      nc - number of components per grid point
136247c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
136347c6ae99SBarry Smith   */
13641321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
136547c6ae99SBarry Smith   col  = 2*s + 1;
136647c6ae99SBarry Smith 
1367aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1368aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
136947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
137047c6ae99SBarry Smith 
1371785e854fSJed Brown   ierr = PetscMalloc1(col*col*nc*nc,&cols);CHKERRQ(ierr);
137247c6ae99SBarry Smith 
13731411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
13741411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
137547c6ae99SBarry Smith 
137647c6ae99SBarry Smith   /* determine the matrix preallocation information */
137747c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
137847c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
13791321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
13801321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
138147c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
13821321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
13831321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
138447c6ae99SBarry Smith       slot   = i - gxs + gnx*(j - gys);
138547c6ae99SBarry Smith 
138647c6ae99SBarry Smith       /* Find block columns in block row */
138747c6ae99SBarry Smith       cnt = 0;
138847c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
138947c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1390aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
139147c6ae99SBarry Smith             cols[cnt++] = slot + ii + gnx*jj;
139247c6ae99SBarry Smith           }
139347c6ae99SBarry Smith         }
139447c6ae99SBarry Smith       }
1395784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
139647c6ae99SBarry Smith     }
139747c6ae99SBarry Smith   }
139847c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
139947c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
140047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
140147c6ae99SBarry Smith 
1402784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1403784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
140447c6ae99SBarry Smith 
140547c6ae99SBarry Smith   /*
140647c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
140747c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
140847c6ae99SBarry Smith     PETSc ordering.
140947c6ae99SBarry Smith   */
1410fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1411*1795a4d1SJed Brown     ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr);
141247c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
14131321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14141321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
141547c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
14161321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14171321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
141847c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
141947c6ae99SBarry Smith         cnt  = 0;
142047c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
142147c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1422aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
142347c6ae99SBarry Smith               cols[cnt++] = slot + ii + gnx*jj;
142447c6ae99SBarry Smith             }
142547c6ae99SBarry Smith           }
142647c6ae99SBarry Smith         }
142747c6ae99SBarry Smith         ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
142847c6ae99SBarry Smith       }
142947c6ae99SBarry Smith     }
143047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
143147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
143247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
143347c6ae99SBarry Smith   }
143447c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
143547c6ae99SBarry Smith   PetscFunctionReturn(0);
143647c6ae99SBarry Smith }
143747c6ae99SBarry Smith 
143847c6ae99SBarry Smith #undef __FUNCT__
1439950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIBAIJ"
1440950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM da,Mat J)
144147c6ae99SBarry Smith {
144247c6ae99SBarry Smith   PetscErrorCode         ierr;
144347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
144447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
144547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
144647c6ae99SBarry Smith   MPI_Comm               comm;
144747c6ae99SBarry Smith   PetscScalar            *values;
14481321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1449aa219208SBarry Smith   DMDAStencilType        st;
145047c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
145147c6ae99SBarry Smith 
145247c6ae99SBarry Smith   PetscFunctionBegin;
145347c6ae99SBarry Smith   /*
145447c6ae99SBarry Smith          nc - number of components per grid point
145547c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
145647c6ae99SBarry Smith 
145747c6ae99SBarry Smith   */
14581321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
145947c6ae99SBarry Smith   col  = 2*s + 1;
146047c6ae99SBarry Smith 
1461aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1462aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
146347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
146447c6ae99SBarry Smith 
1465785e854fSJed Brown   ierr = PetscMalloc1(col*col*col,&cols);CHKERRQ(ierr);
146647c6ae99SBarry Smith 
14671411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
14681411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
146947c6ae99SBarry Smith 
147047c6ae99SBarry Smith   /* determine the matrix preallocation information */
147147c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
147247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
14731321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14741321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
147547c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
14761321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14771321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
147847c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
14791321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
14801321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
148147c6ae99SBarry Smith 
148247c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
148347c6ae99SBarry Smith 
148447c6ae99SBarry Smith         /* Find block columns in block row */
148547c6ae99SBarry Smith         cnt = 0;
148647c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
148747c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
148847c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1489aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
149047c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
149147c6ae99SBarry Smith               }
149247c6ae99SBarry Smith             }
149347c6ae99SBarry Smith           }
149447c6ae99SBarry Smith         }
1495784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
149647c6ae99SBarry Smith       }
149747c6ae99SBarry Smith     }
149847c6ae99SBarry Smith   }
149947c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
150047c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
150147c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
150247c6ae99SBarry Smith 
1503784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1504784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
150547c6ae99SBarry Smith 
150647c6ae99SBarry Smith   /*
150747c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
150847c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
150947c6ae99SBarry Smith     PETSc ordering.
151047c6ae99SBarry Smith   */
1511fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1512*1795a4d1SJed Brown     ierr = PetscCalloc1(col*col*col*nc*nc,&values);CHKERRQ(ierr);
151347c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
15141321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
15151321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
151647c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
15171321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15181321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
151947c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
15201321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15211321219cSEthan Coon           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
152247c6ae99SBarry Smith 
152347c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
152447c6ae99SBarry Smith 
152547c6ae99SBarry Smith           cnt = 0;
152647c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
152747c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
152847c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1529aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
153047c6ae99SBarry Smith                   cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
153147c6ae99SBarry Smith                 }
153247c6ae99SBarry Smith               }
153347c6ae99SBarry Smith             }
153447c6ae99SBarry Smith           }
153547c6ae99SBarry Smith           ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
153647c6ae99SBarry Smith         }
153747c6ae99SBarry Smith       }
153847c6ae99SBarry Smith     }
153947c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
154047c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154147c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
154247c6ae99SBarry Smith   }
154347c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
154447c6ae99SBarry Smith   PetscFunctionReturn(0);
154547c6ae99SBarry Smith }
154647c6ae99SBarry Smith 
154747c6ae99SBarry Smith #undef __FUNCT__
154847c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular"
154947c6ae99SBarry Smith /*
155047c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
155147c6ae99SBarry Smith   identify in the local ordering with periodic domain.
155247c6ae99SBarry Smith */
155347c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
155447c6ae99SBarry Smith {
155547c6ae99SBarry Smith   PetscErrorCode ierr;
155647c6ae99SBarry Smith   PetscInt       i,n;
155747c6ae99SBarry Smith 
155847c6ae99SBarry Smith   PetscFunctionBegin;
155947c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr);
156047c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr);
156147c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
156247c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
156347c6ae99SBarry Smith   }
156447c6ae99SBarry Smith   *cnt = n;
156547c6ae99SBarry Smith   PetscFunctionReturn(0);
156647c6ae99SBarry Smith }
156747c6ae99SBarry Smith 
156847c6ae99SBarry Smith #undef __FUNCT__
1569950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPISBAIJ"
1570950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM da,Mat J)
157147c6ae99SBarry Smith {
157247c6ae99SBarry Smith   PetscErrorCode         ierr;
157347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
157447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
157547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
157647c6ae99SBarry Smith   MPI_Comm               comm;
157747c6ae99SBarry Smith   PetscScalar            *values;
15781321219cSEthan Coon   DMDABoundaryType       bx,by;
1579aa219208SBarry Smith   DMDAStencilType        st;
158047c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
158147c6ae99SBarry Smith 
158247c6ae99SBarry Smith   PetscFunctionBegin;
158347c6ae99SBarry Smith   /*
158447c6ae99SBarry Smith      nc - number of components per grid point
158547c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
158647c6ae99SBarry Smith   */
15871321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
158847c6ae99SBarry Smith   col  = 2*s + 1;
158947c6ae99SBarry Smith 
1590aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1591aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
159247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
159347c6ae99SBarry Smith 
1594785e854fSJed Brown   ierr = PetscMalloc1(col*col*nc*nc,&cols);CHKERRQ(ierr);
159547c6ae99SBarry Smith 
15961411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
15971411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
159847c6ae99SBarry Smith 
159947c6ae99SBarry Smith   /* determine the matrix preallocation information */
1600eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
160147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
16021321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16031321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
160447c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
16051321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16061321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
160747c6ae99SBarry Smith       slot   = i - gxs + gnx*(j - gys);
160847c6ae99SBarry Smith 
160947c6ae99SBarry Smith       /* Find block columns in block row */
161047c6ae99SBarry Smith       cnt = 0;
161147c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
161247c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1613aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) {
161447c6ae99SBarry Smith             cols[cnt++] = slot + ii + gnx*jj;
161547c6ae99SBarry Smith           }
161647c6ae99SBarry Smith         }
161747c6ae99SBarry Smith       }
161847c6ae99SBarry Smith       ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
161947c6ae99SBarry Smith       ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
162047c6ae99SBarry Smith     }
162147c6ae99SBarry Smith   }
162247c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
162347c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
162447c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
162547c6ae99SBarry Smith 
1626784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1627784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
162847c6ae99SBarry Smith 
162947c6ae99SBarry Smith   /*
163047c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
163147c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
163247c6ae99SBarry Smith     PETSc ordering.
163347c6ae99SBarry Smith   */
1634fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1635*1795a4d1SJed Brown     ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr);
163647c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
16371321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16381321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
163947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
16401321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16411321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
164247c6ae99SBarry Smith         slot   = i - gxs + gnx*(j - gys);
164347c6ae99SBarry Smith 
164447c6ae99SBarry Smith         /* Find block columns in block row */
164547c6ae99SBarry Smith         cnt = 0;
164647c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
164747c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1648aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) {
164947c6ae99SBarry Smith               cols[cnt++] = slot + ii + gnx*jj;
165047c6ae99SBarry Smith             }
165147c6ae99SBarry Smith           }
165247c6ae99SBarry Smith         }
165347c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
165447c6ae99SBarry Smith         ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
165547c6ae99SBarry Smith       }
165647c6ae99SBarry Smith     }
165747c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
165847c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
165947c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166047c6ae99SBarry Smith   }
166147c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
166247c6ae99SBarry Smith   PetscFunctionReturn(0);
166347c6ae99SBarry Smith }
166447c6ae99SBarry Smith 
166547c6ae99SBarry Smith #undef __FUNCT__
1666950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPISBAIJ"
1667950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM da,Mat J)
166847c6ae99SBarry Smith {
166947c6ae99SBarry Smith   PetscErrorCode         ierr;
167047c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
167147c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
167247c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
167347c6ae99SBarry Smith   MPI_Comm               comm;
167447c6ae99SBarry Smith   PetscScalar            *values;
16751321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1676aa219208SBarry Smith   DMDAStencilType        st;
167747c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
167847c6ae99SBarry Smith 
167947c6ae99SBarry Smith   PetscFunctionBegin;
168047c6ae99SBarry Smith   /*
168147c6ae99SBarry Smith      nc - number of components per grid point
168247c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
168347c6ae99SBarry Smith   */
16841321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
168547c6ae99SBarry Smith   col  = 2*s + 1;
168647c6ae99SBarry Smith 
1687aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1688aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
168947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
169047c6ae99SBarry Smith 
169147c6ae99SBarry Smith   /* create the matrix */
1692785e854fSJed Brown   ierr = PetscMalloc1(col*col*col,&cols);CHKERRQ(ierr);
169347c6ae99SBarry Smith 
16941411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
16951411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
169647c6ae99SBarry Smith 
169747c6ae99SBarry Smith   /* determine the matrix preallocation information */
1698eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
169947c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
17001321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
17011321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
170247c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
17031321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
17041321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
170547c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
17061321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
17071321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
170847c6ae99SBarry Smith 
170947c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
171047c6ae99SBarry Smith 
171147c6ae99SBarry Smith         /* Find block columns in block row */
171247c6ae99SBarry Smith         cnt = 0;
171347c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
171447c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
171547c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1716aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
171747c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
171847c6ae99SBarry Smith               }
171947c6ae99SBarry Smith             }
172047c6ae99SBarry Smith           }
172147c6ae99SBarry Smith         }
172247c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
172347c6ae99SBarry Smith         ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
172447c6ae99SBarry Smith       }
172547c6ae99SBarry Smith     }
172647c6ae99SBarry Smith   }
172747c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
172847c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
172947c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
173047c6ae99SBarry Smith 
1731784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1732784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
173347c6ae99SBarry Smith 
173447c6ae99SBarry Smith   /*
173547c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
173647c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
173747c6ae99SBarry Smith     PETSc ordering.
173847c6ae99SBarry Smith   */
1739fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1740*1795a4d1SJed Brown     ierr = PetscCalloc1(col*col*col*nc*nc,&values);CHKERRQ(ierr);
174147c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
17421321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
17431321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
174447c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
17451321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
17461321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
174747c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
17481321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
17491321219cSEthan Coon           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
175047c6ae99SBarry Smith 
175147c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
175247c6ae99SBarry Smith 
175347c6ae99SBarry Smith           cnt = 0;
175447c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
175547c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
175647c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1757aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
175847c6ae99SBarry Smith                   cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
175947c6ae99SBarry Smith                 }
176047c6ae99SBarry Smith               }
176147c6ae99SBarry Smith             }
176247c6ae99SBarry Smith           }
176347c6ae99SBarry Smith           ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
176447c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
176547c6ae99SBarry Smith         }
176647c6ae99SBarry Smith       }
176747c6ae99SBarry Smith     }
176847c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
176947c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177047c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177147c6ae99SBarry Smith   }
177247c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
177347c6ae99SBarry Smith   PetscFunctionReturn(0);
177447c6ae99SBarry Smith }
177547c6ae99SBarry Smith 
177647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
177747c6ae99SBarry Smith 
177847c6ae99SBarry Smith #undef __FUNCT__
1779950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ_Fill"
1780950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J)
178147c6ae99SBarry Smith {
178247c6ae99SBarry Smith   PetscErrorCode         ierr;
178347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
178447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz;
178547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
178647c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
178747c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
178847c6ae99SBarry Smith   MPI_Comm               comm;
178947c6ae99SBarry Smith   PetscScalar            *values;
17901321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
179147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1792aa219208SBarry Smith   DMDAStencilType        st;
179347c6ae99SBarry Smith 
179447c6ae99SBarry Smith   PetscFunctionBegin;
179547c6ae99SBarry Smith   /*
179647c6ae99SBarry Smith          nc - number of components per grid point
179747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
179847c6ae99SBarry Smith 
179947c6ae99SBarry Smith   */
18001321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
180147c6ae99SBarry Smith   col  = 2*s + 1;
1802ce94432eSBarry Smith   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
180347c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1804ce94432eSBarry Smith   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
180547c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1806ce94432eSBarry Smith   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
180747c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
180847c6ae99SBarry Smith 
1809aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1810aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
181147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
181247c6ae99SBarry Smith 
1813785e854fSJed Brown   ierr = PetscMalloc1(col*col*col*nc,&cols);CHKERRQ(ierr);
18141411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
18151411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
181647c6ae99SBarry Smith 
181747c6ae99SBarry Smith   /* determine the matrix preallocation information */
181847c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
181947c6ae99SBarry Smith 
182047c6ae99SBarry Smith 
182147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
18221321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
18231321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
182447c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
18251321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
18261321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
182747c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
18281321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
18291321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
183047c6ae99SBarry Smith 
183147c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
183247c6ae99SBarry Smith 
183347c6ae99SBarry Smith         for (l=0; l<nc; l++) {
183447c6ae99SBarry Smith           cnt = 0;
183547c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
183647c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
183747c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
183847c6ae99SBarry Smith                 if (ii || jj || kk) {
1839aa219208SBarry Smith                   if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
18408865f1eaSKarl 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);
184147c6ae99SBarry Smith                   }
184247c6ae99SBarry Smith                 } else {
184347c6ae99SBarry Smith                   if (dfill) {
18448865f1eaSKarl 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);
184547c6ae99SBarry Smith                   } else {
18468865f1eaSKarl Rupp                     for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
184747c6ae99SBarry Smith                   }
184847c6ae99SBarry Smith                 }
184947c6ae99SBarry Smith               }
185047c6ae99SBarry Smith             }
185147c6ae99SBarry Smith           }
185247c6ae99SBarry Smith           row  = l + nc*(slot);
1853784ac674SJed Brown           ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
185447c6ae99SBarry Smith         }
185547c6ae99SBarry Smith       }
185647c6ae99SBarry Smith     }
185747c6ae99SBarry Smith   }
185847c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
185947c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
186047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1861784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1862784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
186347c6ae99SBarry Smith 
186447c6ae99SBarry Smith   /*
186547c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
186647c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
186747c6ae99SBarry Smith     PETSc ordering.
186847c6ae99SBarry Smith   */
1869fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1870*1795a4d1SJed Brown     ierr = PetscCalloc1(col*col*col*nc*nc*nc,&values);CHKERRQ(ierr);
187147c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
18721321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
18731321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
187447c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
18751321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
18761321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
187747c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
18781321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
18791321219cSEthan Coon           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
188047c6ae99SBarry Smith 
188147c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
188247c6ae99SBarry Smith 
188347c6ae99SBarry Smith           for (l=0; l<nc; l++) {
188447c6ae99SBarry Smith             cnt = 0;
188547c6ae99SBarry Smith             for (ii=istart; ii<iend+1; ii++) {
188647c6ae99SBarry Smith               for (jj=jstart; jj<jend+1; jj++) {
188747c6ae99SBarry Smith                 for (kk=kstart; kk<kend+1; kk++) {
188847c6ae99SBarry Smith                   if (ii || jj || kk) {
1889aa219208SBarry Smith                     if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
18908865f1eaSKarl 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);
189147c6ae99SBarry Smith                     }
189247c6ae99SBarry Smith                   } else {
189347c6ae99SBarry Smith                     if (dfill) {
18948865f1eaSKarl 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);
189547c6ae99SBarry Smith                     } else {
18968865f1eaSKarl Rupp                       for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
189747c6ae99SBarry Smith                     }
189847c6ae99SBarry Smith                   }
189947c6ae99SBarry Smith                 }
190047c6ae99SBarry Smith               }
190147c6ae99SBarry Smith             }
190247c6ae99SBarry Smith             row  = l + nc*(slot);
190347c6ae99SBarry Smith             ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
190447c6ae99SBarry Smith           }
190547c6ae99SBarry Smith         }
190647c6ae99SBarry Smith       }
190747c6ae99SBarry Smith     }
190847c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
190947c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
191047c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
191147c6ae99SBarry Smith   }
191247c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
191347c6ae99SBarry Smith   PetscFunctionReturn(0);
191447c6ae99SBarry Smith }
1915