xref: /petsc/src/dm/impls/da/fdda.c (revision b2566f29c2b6470df769aa9f7deb9e2726b0959e)
147c6ae99SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I      "petscdmda.h"     I*/
307475bc1SBarry Smith #include <petscmat.h>
4af0996ceSBarry Smith #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   }
34854ce69bSBarry Smith   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 */
1021795a4d1SJed 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;
121bff4a2f0SMatthew G. Knepley   DMBoundaryType 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) {
158bff4a2f0SMatthew G. Knepley       if ((m==1 && bx == DM_BOUNDARY_PERIODIC) || (n==1 && by == DM_BOUNDARY_PERIODIC) || (p==1 && bz == DM_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;
209bff4a2f0SMatthew G. Knepley   DMBoundaryType 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 
231bff4a2f0SMatthew G. Knepley     if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\
23247c6ae99SBarry Smith                                                             by 2*stencil_width + 1 (%d)\n", m, col);
233bff4a2f0SMatthew G. Knepley     if (by == DM_BOUNDARY_PERIODIC && (n % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\
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));
247aaf3ff59SMatthew G. Knepley         ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,PETSC_OWN_POINTER,&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));
263aaf3ff59SMatthew G. Knepley         ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,PETSC_OWN_POINTER,&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;
285bff4a2f0SMatthew G. Knepley   DMBoundaryType 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;
298bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
29947c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
300bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC && (n % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
30147c6ae99SBarry Smith                                                          by 2*stencil_width + 1\n");
302bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC && (p % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
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));
324aaf3ff59SMatthew G. Knepley       ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,PETSC_OWN_POINTER,&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));
342aaf3ff59SMatthew G. Knepley       ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,PETSC_OWN_POINTER,&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;
361bff4a2f0SMatthew G. Knepley   DMBoundaryType 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 
374bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ2(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points %d is divisible\n\
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       }
408aaf3ff59SMatthew G. Knepley       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,PETSC_OWN_POINTER,&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);
422aaf3ff59SMatthew G. Knepley       ierr    = ISColoringCreate(comm,ncolors,nc*gnx,colors,PETSC_OWN_POINTER,&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;
439bff4a2f0SMatthew G. Knepley   DMBoundaryType 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 
454bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && (m % 5)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible by 5\n");
455bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC && (n % 5)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible by 5\n");
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;
470aaf3ff59SMatthew G. Knepley       ierr    = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,PETSC_OWN_POINTER,&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;
485aaf3ff59SMatthew G. Knepley       ierr    = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,PETSC_OWN_POINTER,&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);
555854ce69bSBarry Smith   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);
566162ae106SBarry Smith   ierr = (Anatural->ops->view)(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);
597854ce69bSBarry Smith   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);
60428be2f97SBarry 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;
673*b2566f29SBarry Smith         ierr    = MPIU_Allreduce(&bsLocal, &bs, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)da));CHKERRQ(ierr);
67437d0c07bSMatthew G Knepley       }
6751795a4d1SJed 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;
705c73cfb54SMatthew G. Knepley   dim = da->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 {
77545b6f7e9SBarry Smith     ISLocalToGlobalMapping ltog;
776869776cdSLisandro Dalcin     ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
7772949035bSJed Brown     ierr = MatSetUp(A);CHKERRQ(ierr);
778869776cdSLisandro Dalcin     ierr = MatSetLocalToGlobalMapping(A,ltog,ltog);CHKERRQ(ierr);
77947c6ae99SBarry Smith   }
780aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
78147c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
782c688c046SMatthew G Knepley   ierr = MatSetDM(A,da);CHKERRQ(ierr);
78347c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
78447c6ae99SBarry Smith   if (size > 1) {
78547c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
78647c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void))MatView_MPI_DA);CHKERRQ(ierr);
78747c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void))MatLoad_MPI_DA);CHKERRQ(ierr);
78847c6ae99SBarry Smith   }
78947c6ae99SBarry Smith   *J = A;
79047c6ae99SBarry Smith   PetscFunctionReturn(0);
79147c6ae99SBarry Smith }
79247c6ae99SBarry Smith 
79347c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
79447c6ae99SBarry Smith #undef __FUNCT__
795950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ"
796950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ(DM da,Mat J)
79747c6ae99SBarry Smith {
79847c6ae99SBarry Smith   PetscErrorCode         ierr;
7990298fd71SBarry 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;
80047c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
80147c6ae99SBarry Smith   MPI_Comm               comm;
80247c6ae99SBarry Smith   PetscScalar            *values;
803bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
80445b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
805aa219208SBarry Smith   DMDAStencilType        st;
80647c6ae99SBarry Smith 
80747c6ae99SBarry Smith   PetscFunctionBegin;
80847c6ae99SBarry Smith   /*
80947c6ae99SBarry Smith          nc - number of components per grid point
81047c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
81147c6ae99SBarry Smith 
81247c6ae99SBarry Smith   */
8131321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
81447c6ae99SBarry Smith   col  = 2*s + 1;
815aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
816aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
81747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
81847c6ae99SBarry Smith 
819dcca6d9dSJed Brown   ierr = PetscMalloc2(nc,&rows,col*col*nc*nc,&cols);CHKERRQ(ierr);
8201411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
82147c6ae99SBarry Smith 
82206ca8cadSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
82347c6ae99SBarry Smith   /* determine the matrix preallocation information */
82447c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
82547c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
82647c6ae99SBarry Smith 
827bff4a2f0SMatthew G. Knepley     pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
828bff4a2f0SMatthew G. Knepley     pend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
82947c6ae99SBarry Smith 
83047c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
83147c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
83247c6ae99SBarry Smith 
833bff4a2f0SMatthew G. Knepley       lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
834bff4a2f0SMatthew G. Knepley       lend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
83547c6ae99SBarry Smith 
83647c6ae99SBarry Smith       cnt = 0;
83747c6ae99SBarry Smith       for (k=0; k<nc; k++) {
83847c6ae99SBarry Smith         for (l=lstart; l<lend+1; l++) {
83947c6ae99SBarry Smith           for (p=pstart; p<pend+1; p++) {
840aa219208SBarry Smith             if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
84147c6ae99SBarry Smith               cols[cnt++] = k + nc*(slot + gnx*l + p);
84247c6ae99SBarry Smith             }
84347c6ae99SBarry Smith           }
84447c6ae99SBarry Smith         }
84547c6ae99SBarry Smith         rows[k] = k + nc*(slot);
84647c6ae99SBarry Smith       }
847784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
84847c6ae99SBarry Smith     }
84947c6ae99SBarry Smith   }
850f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
85147c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
85247c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
85347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
85447c6ae99SBarry Smith 
855784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
85647c6ae99SBarry Smith 
85747c6ae99SBarry Smith   /*
85847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
85947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
86047c6ae99SBarry Smith     PETSc ordering.
86147c6ae99SBarry Smith   */
862fcfd50ebSBarry Smith   if (!da->prealloc_only) {
8631795a4d1SJed Brown     ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr);
86447c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
86547c6ae99SBarry Smith 
866bff4a2f0SMatthew G. Knepley       pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
867bff4a2f0SMatthew G. Knepley       pend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
86847c6ae99SBarry Smith 
86947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
87047c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
87147c6ae99SBarry Smith 
872bff4a2f0SMatthew G. Knepley         lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
873bff4a2f0SMatthew G. Knepley         lend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
87447c6ae99SBarry Smith 
87547c6ae99SBarry Smith         cnt = 0;
87647c6ae99SBarry Smith         for (k=0; k<nc; k++) {
87747c6ae99SBarry Smith           for (l=lstart; l<lend+1; l++) {
87847c6ae99SBarry Smith             for (p=pstart; p<pend+1; p++) {
879aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
88047c6ae99SBarry Smith                 cols[cnt++] = k + nc*(slot + gnx*l + p);
88147c6ae99SBarry Smith               }
88247c6ae99SBarry Smith             }
88347c6ae99SBarry Smith           }
88447c6ae99SBarry Smith           rows[k] = k + nc*(slot);
88547c6ae99SBarry Smith         }
88647c6ae99SBarry Smith         ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
88747c6ae99SBarry Smith       }
88847c6ae99SBarry Smith     }
88947c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
89047c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89147c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
892189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
89347c6ae99SBarry Smith   }
89447c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
89547c6ae99SBarry Smith   PetscFunctionReturn(0);
89647c6ae99SBarry Smith }
89747c6ae99SBarry Smith 
89847c6ae99SBarry Smith #undef __FUNCT__
899950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIAIJ_Fill"
900950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J)
90147c6ae99SBarry Smith {
90247c6ae99SBarry Smith   PetscErrorCode         ierr;
90347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
904c0ab637bSBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,maxcnt = 0,l,p;
90547c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
90647c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
90747c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
90847c6ae99SBarry Smith   MPI_Comm               comm;
90947c6ae99SBarry Smith   PetscScalar            *values;
910bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
91145b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
912aa219208SBarry Smith   DMDAStencilType        st;
91347c6ae99SBarry Smith 
91447c6ae99SBarry Smith   PetscFunctionBegin;
91547c6ae99SBarry Smith   /*
91647c6ae99SBarry Smith          nc - number of components per grid point
91747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
91847c6ae99SBarry Smith 
91947c6ae99SBarry Smith   */
9201321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
92147c6ae99SBarry Smith   col  = 2*s + 1;
922aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
923aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
92447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
92547c6ae99SBarry Smith 
9264b26d1cfSBarry Smith   ierr = PetscMalloc1(col*col*nc,&cols);CHKERRQ(ierr);
9271411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
92847c6ae99SBarry Smith 
92906ca8cadSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
93047c6ae99SBarry Smith   /* determine the matrix preallocation information */
93147c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
93247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
93347c6ae99SBarry Smith 
934bff4a2f0SMatthew G. Knepley     pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
935bff4a2f0SMatthew G. Knepley     pend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
93647c6ae99SBarry Smith 
93747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
93847c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
93947c6ae99SBarry Smith 
940bff4a2f0SMatthew G. Knepley       lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
941bff4a2f0SMatthew G. Knepley       lend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
94247c6ae99SBarry Smith 
94347c6ae99SBarry Smith       for (k=0; k<nc; k++) {
94447c6ae99SBarry Smith         cnt = 0;
94547c6ae99SBarry Smith         for (l=lstart; l<lend+1; l++) {
94647c6ae99SBarry Smith           for (p=pstart; p<pend+1; p++) {
94747c6ae99SBarry Smith             if (l || p) {
948aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
9498865f1eaSKarl Rupp                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p);
95047c6ae99SBarry Smith               }
95147c6ae99SBarry Smith             } else {
95247c6ae99SBarry Smith               if (dfill) {
9538865f1eaSKarl Rupp                 for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p);
95447c6ae99SBarry Smith               } else {
9558865f1eaSKarl Rupp                 for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + gnx*l + p);
95647c6ae99SBarry Smith               }
95747c6ae99SBarry Smith             }
95847c6ae99SBarry Smith           }
95947c6ae99SBarry Smith         }
96047c6ae99SBarry Smith         row    = k + nc*(slot);
961c0ab637bSBarry Smith         maxcnt = PetscMax(maxcnt,cnt);
962784ac674SJed Brown         ierr   = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
96347c6ae99SBarry Smith       }
96447c6ae99SBarry Smith     }
96547c6ae99SBarry Smith   }
96647c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
96747c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
96847c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
969784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
97047c6ae99SBarry Smith 
97147c6ae99SBarry Smith   /*
97247c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
97347c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
97447c6ae99SBarry Smith     PETSc ordering.
97547c6ae99SBarry Smith   */
976fcfd50ebSBarry Smith   if (!da->prealloc_only) {
977c0ab637bSBarry Smith     ierr = PetscCalloc1(maxcnt,&values);CHKERRQ(ierr);
97847c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
97947c6ae99SBarry Smith 
980bff4a2f0SMatthew G. Knepley       pstart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
981bff4a2f0SMatthew G. Knepley       pend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
98247c6ae99SBarry Smith 
98347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
98447c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
98547c6ae99SBarry Smith 
986bff4a2f0SMatthew G. Knepley         lstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
987bff4a2f0SMatthew G. Knepley         lend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
98847c6ae99SBarry Smith 
98947c6ae99SBarry Smith         for (k=0; k<nc; k++) {
99047c6ae99SBarry Smith           cnt = 0;
99147c6ae99SBarry Smith           for (l=lstart; l<lend+1; l++) {
99247c6ae99SBarry Smith             for (p=pstart; p<pend+1; p++) {
99347c6ae99SBarry Smith               if (l || p) {
994aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
9958865f1eaSKarl Rupp                   for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++) cols[cnt++] = ofill[ifill_col] + nc*(slot + gnx*l + p);
99647c6ae99SBarry Smith                 }
99747c6ae99SBarry Smith               } else {
99847c6ae99SBarry Smith                 if (dfill) {
9998865f1eaSKarl Rupp                   for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++) cols[cnt++] = dfill[ifill_col] + nc*(slot + gnx*l + p);
100047c6ae99SBarry Smith                 } else {
10018865f1eaSKarl Rupp                   for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + gnx*l + p);
100247c6ae99SBarry Smith                 }
100347c6ae99SBarry Smith               }
100447c6ae99SBarry Smith             }
100547c6ae99SBarry Smith           }
100647c6ae99SBarry Smith           row  = k + nc*(slot);
100747c6ae99SBarry Smith           ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
100847c6ae99SBarry Smith         }
100947c6ae99SBarry Smith       }
101047c6ae99SBarry Smith     }
101147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
101247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
101347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1014189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
101547c6ae99SBarry Smith   }
101647c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
101747c6ae99SBarry Smith   PetscFunctionReturn(0);
101847c6ae99SBarry Smith }
101947c6ae99SBarry Smith 
102047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
102147c6ae99SBarry Smith 
102247c6ae99SBarry Smith #undef __FUNCT__
1023950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ"
1024950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ(DM da,Mat J)
102547c6ae99SBarry Smith {
102647c6ae99SBarry Smith   PetscErrorCode         ierr;
102747c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
10280298fd71SBarry Smith   PetscInt               m,n,dim,s,*cols = NULL,k,nc,*rows = NULL,col,cnt,l,p,*dnz = NULL,*onz = NULL;
102947c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
103047c6ae99SBarry Smith   MPI_Comm               comm;
103147c6ae99SBarry Smith   PetscScalar            *values;
1032bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
103345b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
1034aa219208SBarry Smith   DMDAStencilType        st;
103547c6ae99SBarry Smith 
103647c6ae99SBarry Smith   PetscFunctionBegin;
103747c6ae99SBarry Smith   /*
103847c6ae99SBarry Smith          nc - number of components per grid point
103947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
104047c6ae99SBarry Smith 
104147c6ae99SBarry Smith   */
10421321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
104347c6ae99SBarry Smith   col  = 2*s + 1;
104447c6ae99SBarry Smith 
1045aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1046aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
104747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
104847c6ae99SBarry Smith 
1049dcca6d9dSJed Brown   ierr = PetscMalloc2(nc,&rows,col*col*col*nc*nc,&cols);CHKERRQ(ierr);
10501411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
105147c6ae99SBarry Smith 
105206ca8cadSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
105347c6ae99SBarry Smith   /* determine the matrix preallocation information */
105447c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
105547c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1056bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1057bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
105847c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1059bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1060bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
106147c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1062bff4a2f0SMatthew G. Knepley         kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1063bff4a2f0SMatthew G. Knepley         kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
106447c6ae99SBarry Smith 
106547c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
106647c6ae99SBarry Smith 
106747c6ae99SBarry Smith         cnt = 0;
106847c6ae99SBarry Smith         for (l=0; l<nc; l++) {
106947c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
107047c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
107147c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1072aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
107347c6ae99SBarry Smith                   cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
107447c6ae99SBarry Smith                 }
107547c6ae99SBarry Smith               }
107647c6ae99SBarry Smith             }
107747c6ae99SBarry Smith           }
107847c6ae99SBarry Smith           rows[l] = l + nc*(slot);
107947c6ae99SBarry Smith         }
1080784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
108147c6ae99SBarry Smith       }
108247c6ae99SBarry Smith     }
108347c6ae99SBarry Smith   }
1084f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
108547c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
108647c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
108747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1088784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
108947c6ae99SBarry Smith 
109047c6ae99SBarry Smith   /*
109147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
109247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
109347c6ae99SBarry Smith     PETSc ordering.
109447c6ae99SBarry Smith   */
1095fcfd50ebSBarry Smith   if (!da->prealloc_only) {
10961795a4d1SJed Brown     ierr = PetscCalloc1(col*col*col*nc*nc*nc,&values);CHKERRQ(ierr);
109747c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1098bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1099bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
110047c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1101bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1102bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
110347c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
1104bff4a2f0SMatthew G. Knepley           kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1105bff4a2f0SMatthew G. Knepley           kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
110647c6ae99SBarry Smith 
110747c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
110847c6ae99SBarry Smith 
110947c6ae99SBarry Smith           cnt = 0;
111047c6ae99SBarry Smith           for (l=0; l<nc; l++) {
111147c6ae99SBarry Smith             for (ii=istart; ii<iend+1; ii++) {
111247c6ae99SBarry Smith               for (jj=jstart; jj<jend+1; jj++) {
111347c6ae99SBarry Smith                 for (kk=kstart; kk<kend+1; kk++) {
1114aa219208SBarry Smith                   if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
111547c6ae99SBarry Smith                     cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
111647c6ae99SBarry Smith                   }
111747c6ae99SBarry Smith                 }
111847c6ae99SBarry Smith               }
111947c6ae99SBarry Smith             }
112047c6ae99SBarry Smith             rows[l] = l + nc*(slot);
112147c6ae99SBarry Smith           }
112247c6ae99SBarry Smith           ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
112347c6ae99SBarry Smith         }
112447c6ae99SBarry Smith       }
112547c6ae99SBarry Smith     }
112647c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
112747c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
112847c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1129189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
113047c6ae99SBarry Smith   }
113147c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
113247c6ae99SBarry Smith   PetscFunctionReturn(0);
113347c6ae99SBarry Smith }
113447c6ae99SBarry Smith 
113547c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
113647c6ae99SBarry Smith 
113747c6ae99SBarry Smith #undef __FUNCT__
1138ce308e1dSBarry Smith #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ_Fill"
1139ce308e1dSBarry Smith PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ_Fill(DM da,Mat J)
1140ce308e1dSBarry Smith {
1141ce308e1dSBarry Smith   PetscErrorCode         ierr;
1142ce308e1dSBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1143ce308e1dSBarry Smith   PetscInt               xs,nx,i,j,gxs,gnx,row,k,l;
11448d4c968fSBarry Smith   PetscInt               m,dim,s,*cols = NULL,nc,cnt,maxcnt = 0,*ocols;
11450acb5bebSBarry Smith   PetscInt               *ofill = dd->ofill,*dfill = dd->dfill;
1146ce308e1dSBarry Smith   PetscScalar            *values;
1147bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
114845b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
1149ce308e1dSBarry Smith   PetscMPIInt            rank,size;
1150ce308e1dSBarry Smith 
1151ce308e1dSBarry Smith   PetscFunctionBegin;
1152bff4a2f0SMatthew G. Knepley   if (dd->bx == DM_BOUNDARY_PERIODIC) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"With fill provided not implemented with periodic boundary conditions");
1153ce94432eSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)da),&rank);CHKERRQ(ierr);
1154ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)da),&size);CHKERRQ(ierr);
1155ce308e1dSBarry Smith 
1156ce308e1dSBarry Smith   /*
1157ce308e1dSBarry Smith          nc - number of components per grid point
1158ce308e1dSBarry Smith 
1159ce308e1dSBarry Smith   */
1160ce308e1dSBarry Smith   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
1161ce308e1dSBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1162ce308e1dSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
1163ce308e1dSBarry Smith 
1164ce308e1dSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
11651795a4d1SJed Brown   ierr = PetscCalloc2(nx*nc,&cols,nx*nc,&ocols);CHKERRQ(ierr);
1166ce308e1dSBarry Smith 
1167554c65c0SBarry Smith   if (nx < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Need at least two grid points per process");
1168ce308e1dSBarry Smith   /*
1169ce308e1dSBarry Smith         note should be smaller for first and last process with no periodic
1170ce308e1dSBarry Smith         does not handle dfill
1171ce308e1dSBarry Smith   */
1172ce308e1dSBarry Smith   cnt = 0;
1173ce308e1dSBarry Smith   /* coupling with process to the left */
1174ce308e1dSBarry Smith   for (i=0; i<s; i++) {
1175ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
1176ce308e1dSBarry Smith       ocols[cnt] = ((!rank) ? 0 : (s - i)*(ofill[j+1] - ofill[j]));
11770acb5bebSBarry Smith       cols[cnt]  = dfill[j+1] - dfill[j] + (s + i)*(ofill[j+1] - ofill[j]);
1178c0ab637bSBarry Smith       maxcnt = PetscMax(maxcnt,ocols[cnt]+cols[cnt]);
1179ce308e1dSBarry Smith       cnt++;
1180ce308e1dSBarry Smith     }
1181ce308e1dSBarry Smith   }
1182ce308e1dSBarry Smith   for (i=s; i<nx-s; i++) {
1183ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
11840acb5bebSBarry Smith       cols[cnt] = dfill[j+1] - dfill[j] + 2*s*(ofill[j+1] - ofill[j]);
1185c0ab637bSBarry Smith       maxcnt = PetscMax(maxcnt,ocols[cnt]+cols[cnt]);
1186ce308e1dSBarry Smith       cnt++;
1187ce308e1dSBarry Smith     }
1188ce308e1dSBarry Smith   }
1189ce308e1dSBarry Smith   /* coupling with process to the right */
1190ce308e1dSBarry Smith   for (i=nx-s; i<nx; i++) {
1191ce308e1dSBarry Smith     for (j=0; j<nc; j++) {
1192ce308e1dSBarry Smith       ocols[cnt] = ((rank == (size-1)) ? 0 : (i - nx + s + 1)*(ofill[j+1] - ofill[j]));
11930acb5bebSBarry Smith       cols[cnt]  = dfill[j+1] - dfill[j] + (s + nx - i - 1)*(ofill[j+1] - ofill[j]);
1194c0ab637bSBarry Smith       maxcnt = PetscMax(maxcnt,ocols[cnt]+cols[cnt]);
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 = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1205ce308e1dSBarry Smith 
1206ce308e1dSBarry Smith   /*
1207ce308e1dSBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1208ce308e1dSBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1209ce308e1dSBarry Smith     PETSc ordering.
1210ce308e1dSBarry Smith   */
1211ce308e1dSBarry Smith   if (!da->prealloc_only) {
1212c0ab637bSBarry Smith     ierr = PetscCalloc2(maxcnt,&values,maxcnt,&cols);CHKERRQ(ierr);
1213ce308e1dSBarry Smith 
1214ce308e1dSBarry Smith     row = xs*nc;
1215ce308e1dSBarry Smith     /* coupling with process to the left */
1216ce308e1dSBarry Smith     for (i=xs; i<xs+s; i++) {
1217ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1218ce308e1dSBarry Smith         cnt = 0;
1219ce308e1dSBarry Smith         if (rank) {
1220ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1221ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1222ce308e1dSBarry Smith           }
1223ce308e1dSBarry Smith         }
12240acb5bebSBarry Smith         if (dfill) {
12250acb5bebSBarry Smith           for (k=dfill[j]; k<dfill[j+1]; k++) {
12260acb5bebSBarry Smith             cols[cnt++] = i*nc + dfill[k];
12270acb5bebSBarry Smith           }
12280acb5bebSBarry Smith         } else {
1229ce308e1dSBarry Smith           for (k=0; k<nc; k++) {
1230ce308e1dSBarry Smith             cols[cnt++] = i*nc + k;
1231ce308e1dSBarry Smith           }
12320acb5bebSBarry Smith         }
1233ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1234ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1235ce308e1dSBarry Smith         }
1236ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1237ce308e1dSBarry Smith         row++;
1238ce308e1dSBarry Smith       }
1239ce308e1dSBarry Smith     }
1240ce308e1dSBarry Smith     for (i=xs+s; i<xs+nx-s; i++) {
1241ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1242ce308e1dSBarry Smith         cnt = 0;
1243ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1244ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1245ce308e1dSBarry Smith         }
12460acb5bebSBarry Smith         if (dfill) {
12470acb5bebSBarry Smith           for (k=dfill[j]; k<dfill[j+1]; k++) {
12480acb5bebSBarry Smith             cols[cnt++] = i*nc + dfill[k];
12490acb5bebSBarry Smith           }
12500acb5bebSBarry Smith         } else {
1251ce308e1dSBarry Smith           for (k=0; k<nc; k++) {
1252ce308e1dSBarry Smith             cols[cnt++] = i*nc + k;
1253ce308e1dSBarry Smith           }
12540acb5bebSBarry Smith         }
1255ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1256ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1257ce308e1dSBarry Smith         }
1258ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1259ce308e1dSBarry Smith         row++;
1260ce308e1dSBarry Smith       }
1261ce308e1dSBarry Smith     }
1262ce308e1dSBarry Smith     /* coupling with process to the right */
1263ce308e1dSBarry Smith     for (i=xs+nx-s; i<xs+nx; i++) {
1264ce308e1dSBarry Smith       for (j=0; j<nc; j++) {
1265ce308e1dSBarry Smith         cnt = 0;
1266ce308e1dSBarry Smith         for (l=0; l<s; l++) {
1267ce308e1dSBarry Smith           for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i - s + l)*nc + ofill[k];
1268ce308e1dSBarry Smith         }
12690acb5bebSBarry Smith         if (dfill) {
12700acb5bebSBarry Smith           for (k=dfill[j]; k<dfill[j+1]; k++) {
12710acb5bebSBarry Smith             cols[cnt++] = i*nc + dfill[k];
12720acb5bebSBarry Smith           }
12730acb5bebSBarry Smith         } else {
1274ce308e1dSBarry Smith           for (k=0; k<nc; k++) {
1275ce308e1dSBarry Smith             cols[cnt++] = i*nc + k;
1276ce308e1dSBarry Smith           }
12770acb5bebSBarry Smith         }
1278ce308e1dSBarry Smith         if (rank < size-1) {
1279ce308e1dSBarry Smith           for (l=0; l<s; l++) {
1280ce308e1dSBarry Smith             for (k=ofill[j]; k<ofill[j+1]; k++) cols[cnt++] = (i + s - l)*nc + ofill[k];
1281ce308e1dSBarry Smith           }
1282ce308e1dSBarry Smith         }
1283ce308e1dSBarry Smith         ierr = MatSetValues(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1284ce308e1dSBarry Smith         row++;
1285ce308e1dSBarry Smith       }
1286ce308e1dSBarry Smith     }
1287c0ab637bSBarry Smith     ierr = PetscFree2(values,cols);CHKERRQ(ierr);
1288ce308e1dSBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1289ce308e1dSBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1290189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1291ce308e1dSBarry Smith   }
1292ce308e1dSBarry Smith   PetscFunctionReturn(0);
1293ce308e1dSBarry Smith }
1294ce308e1dSBarry Smith 
1295ce308e1dSBarry Smith /* ---------------------------------------------------------------------------------*/
1296ce308e1dSBarry Smith 
1297ce308e1dSBarry Smith #undef __FUNCT__
1298950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_1d_MPIAIJ"
1299950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_1d_MPIAIJ(DM da,Mat J)
130047c6ae99SBarry Smith {
130147c6ae99SBarry Smith   PetscErrorCode         ierr;
130247c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
13030298fd71SBarry Smith   PetscInt               m,dim,s,*cols = NULL,nc,*rows = NULL,col,cnt,l;
130447c6ae99SBarry Smith   PetscInt               istart,iend;
130547c6ae99SBarry Smith   PetscScalar            *values;
1306bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx;
130745b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
130847c6ae99SBarry Smith 
130947c6ae99SBarry Smith   PetscFunctionBegin;
131047c6ae99SBarry Smith   /*
131147c6ae99SBarry Smith          nc - number of components per grid point
131247c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
131347c6ae99SBarry Smith 
131447c6ae99SBarry Smith   */
13151321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
131647c6ae99SBarry Smith   col  = 2*s + 1;
131747c6ae99SBarry Smith 
1318aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1319aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
132047c6ae99SBarry Smith 
1321f73d5cc4SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
132247c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
132347c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
132447c6ae99SBarry Smith 
13251411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1326784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
132747c6ae99SBarry Smith 
132847c6ae99SBarry Smith   /*
132947c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
133047c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
133147c6ae99SBarry Smith     PETSc ordering.
133247c6ae99SBarry Smith   */
1333fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1334dcca6d9dSJed Brown     ierr = PetscMalloc2(nc,&rows,col*nc*nc,&cols);CHKERRQ(ierr);
13351795a4d1SJed Brown     ierr = PetscCalloc1(col*nc*nc,&values);CHKERRQ(ierr);
133647c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
133747c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
133847c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
133947c6ae99SBarry Smith       slot   = i - gxs;
134047c6ae99SBarry Smith 
134147c6ae99SBarry Smith       cnt = 0;
134247c6ae99SBarry Smith       for (l=0; l<nc; l++) {
134347c6ae99SBarry Smith         for (i1=istart; i1<iend+1; i1++) {
134447c6ae99SBarry Smith           cols[cnt++] = l + nc*(slot + i1);
134547c6ae99SBarry Smith         }
134647c6ae99SBarry Smith         rows[l] = l + nc*(slot);
134747c6ae99SBarry Smith       }
134847c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
134947c6ae99SBarry Smith     }
135047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
135147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
135247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1353189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
135447c6ae99SBarry Smith     ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
1355ce308e1dSBarry Smith   }
135647c6ae99SBarry Smith   PetscFunctionReturn(0);
135747c6ae99SBarry Smith }
135847c6ae99SBarry Smith 
135947c6ae99SBarry Smith #undef __FUNCT__
1360950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPIBAIJ"
1361950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPIBAIJ(DM da,Mat J)
136247c6ae99SBarry Smith {
136347c6ae99SBarry Smith   PetscErrorCode         ierr;
136447c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
136547c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
136647c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
136747c6ae99SBarry Smith   MPI_Comm               comm;
136847c6ae99SBarry Smith   PetscScalar            *values;
1369bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
1370aa219208SBarry Smith   DMDAStencilType        st;
137145b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
137247c6ae99SBarry Smith 
137347c6ae99SBarry Smith   PetscFunctionBegin;
137447c6ae99SBarry Smith   /*
137547c6ae99SBarry Smith      nc - number of components per grid point
137647c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
137747c6ae99SBarry Smith   */
13781321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
137947c6ae99SBarry Smith   col  = 2*s + 1;
138047c6ae99SBarry Smith 
1381aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1382aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
138347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
138447c6ae99SBarry Smith 
1385785e854fSJed Brown   ierr = PetscMalloc1(col*col*nc*nc,&cols);CHKERRQ(ierr);
138647c6ae99SBarry Smith 
13871411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
138847c6ae99SBarry Smith 
138947c6ae99SBarry Smith   /* determine the matrix preallocation information */
139047c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
139147c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1392bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1393bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
139447c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1395bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1396bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
139747c6ae99SBarry Smith       slot   = i - gxs + gnx*(j - gys);
139847c6ae99SBarry Smith 
139947c6ae99SBarry Smith       /* Find block columns in block row */
140047c6ae99SBarry Smith       cnt = 0;
140147c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
140247c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1403aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
140447c6ae99SBarry Smith             cols[cnt++] = slot + ii + gnx*jj;
140547c6ae99SBarry Smith           }
140647c6ae99SBarry Smith         }
140747c6ae99SBarry Smith       }
1408d6e23781SBarry Smith       ierr = MatPreallocateSetLocalBlock(ltog,1,&slot,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
140947c6ae99SBarry Smith     }
141047c6ae99SBarry Smith   }
141147c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
141247c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
141347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
141447c6ae99SBarry Smith 
1415784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
141647c6ae99SBarry Smith 
141747c6ae99SBarry Smith   /*
141847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
141947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
142047c6ae99SBarry Smith     PETSc ordering.
142147c6ae99SBarry Smith   */
1422fcfd50ebSBarry Smith   if (!da->prealloc_only) {
14231795a4d1SJed Brown     ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr);
142447c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1425bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1426bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
142747c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1428bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1429bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
143047c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
143147c6ae99SBarry Smith         cnt  = 0;
143247c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
143347c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1434aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
143547c6ae99SBarry Smith               cols[cnt++] = slot + ii + gnx*jj;
143647c6ae99SBarry Smith             }
143747c6ae99SBarry Smith           }
143847c6ae99SBarry Smith         }
143947c6ae99SBarry Smith         ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
144047c6ae99SBarry Smith       }
144147c6ae99SBarry Smith     }
144247c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
144347c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
144447c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1445189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
144647c6ae99SBarry Smith   }
144747c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
144847c6ae99SBarry Smith   PetscFunctionReturn(0);
144947c6ae99SBarry Smith }
145047c6ae99SBarry Smith 
145147c6ae99SBarry Smith #undef __FUNCT__
1452950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIBAIJ"
1453950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIBAIJ(DM da,Mat J)
145447c6ae99SBarry Smith {
145547c6ae99SBarry Smith   PetscErrorCode         ierr;
145647c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
145747c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
145847c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
145947c6ae99SBarry Smith   MPI_Comm               comm;
146047c6ae99SBarry Smith   PetscScalar            *values;
1461bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
1462aa219208SBarry Smith   DMDAStencilType        st;
146345b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
146447c6ae99SBarry Smith 
146547c6ae99SBarry Smith   PetscFunctionBegin;
146647c6ae99SBarry Smith   /*
146747c6ae99SBarry Smith          nc - number of components per grid point
146847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
146947c6ae99SBarry Smith 
147047c6ae99SBarry Smith   */
14711321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
147247c6ae99SBarry Smith   col  = 2*s + 1;
147347c6ae99SBarry Smith 
1474aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1475aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
147647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
147747c6ae99SBarry Smith 
1478785e854fSJed Brown   ierr = PetscMalloc1(col*col*col,&cols);CHKERRQ(ierr);
147947c6ae99SBarry Smith 
14801411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
148147c6ae99SBarry Smith 
148247c6ae99SBarry Smith   /* determine the matrix preallocation information */
148347c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
148447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1485bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1486bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
148747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1488bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1489bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
149047c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1491bff4a2f0SMatthew G. Knepley         kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1492bff4a2f0SMatthew G. Knepley         kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
149347c6ae99SBarry Smith 
149447c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
149547c6ae99SBarry Smith 
149647c6ae99SBarry Smith         /* Find block columns in block row */
149747c6ae99SBarry Smith         cnt = 0;
149847c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
149947c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
150047c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1501aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
150247c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
150347c6ae99SBarry Smith               }
150447c6ae99SBarry Smith             }
150547c6ae99SBarry Smith           }
150647c6ae99SBarry Smith         }
1507d6e23781SBarry Smith         ierr = MatPreallocateSetLocalBlock(ltog,1,&slot,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
150847c6ae99SBarry Smith       }
150947c6ae99SBarry Smith     }
151047c6ae99SBarry Smith   }
151147c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
151247c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
151347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
151447c6ae99SBarry Smith 
1515784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
151647c6ae99SBarry Smith 
151747c6ae99SBarry Smith   /*
151847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
151947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
152047c6ae99SBarry Smith     PETSc ordering.
152147c6ae99SBarry Smith   */
1522fcfd50ebSBarry Smith   if (!da->prealloc_only) {
15231795a4d1SJed Brown     ierr = PetscCalloc1(col*col*col*nc*nc,&values);CHKERRQ(ierr);
152447c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1525bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1526bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
152747c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1528bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1529bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
153047c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
1531bff4a2f0SMatthew G. Knepley           kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1532bff4a2f0SMatthew G. Knepley           kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
153347c6ae99SBarry Smith 
153447c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
153547c6ae99SBarry Smith 
153647c6ae99SBarry Smith           cnt = 0;
153747c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
153847c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
153947c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1540aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
154147c6ae99SBarry Smith                   cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
154247c6ae99SBarry Smith                 }
154347c6ae99SBarry Smith               }
154447c6ae99SBarry Smith             }
154547c6ae99SBarry Smith           }
154647c6ae99SBarry Smith           ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
154747c6ae99SBarry Smith         }
154847c6ae99SBarry Smith       }
154947c6ae99SBarry Smith     }
155047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
155147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
155247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1553189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
155447c6ae99SBarry Smith   }
155547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
155647c6ae99SBarry Smith   PetscFunctionReturn(0);
155747c6ae99SBarry Smith }
155847c6ae99SBarry Smith 
155947c6ae99SBarry Smith #undef __FUNCT__
156047c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular"
156147c6ae99SBarry Smith /*
156247c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
156347c6ae99SBarry Smith   identify in the local ordering with periodic domain.
156447c6ae99SBarry Smith */
156547c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
156647c6ae99SBarry Smith {
156747c6ae99SBarry Smith   PetscErrorCode ierr;
156847c6ae99SBarry Smith   PetscInt       i,n;
156947c6ae99SBarry Smith 
157047c6ae99SBarry Smith   PetscFunctionBegin;
1571d6e23781SBarry Smith   ierr = ISLocalToGlobalMappingApplyBlock(ltog,1,row,row);CHKERRQ(ierr);
1572d6e23781SBarry Smith   ierr = ISLocalToGlobalMappingApplyBlock(ltog,*cnt,col,col);CHKERRQ(ierr);
157347c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
157447c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
157547c6ae99SBarry Smith   }
157647c6ae99SBarry Smith   *cnt = n;
157747c6ae99SBarry Smith   PetscFunctionReturn(0);
157847c6ae99SBarry Smith }
157947c6ae99SBarry Smith 
158047c6ae99SBarry Smith #undef __FUNCT__
1581950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_2d_MPISBAIJ"
1582950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_2d_MPISBAIJ(DM da,Mat J)
158347c6ae99SBarry Smith {
158447c6ae99SBarry Smith   PetscErrorCode         ierr;
158547c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
158647c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
158747c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
158847c6ae99SBarry Smith   MPI_Comm               comm;
158947c6ae99SBarry Smith   PetscScalar            *values;
1590bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by;
1591aa219208SBarry Smith   DMDAStencilType        st;
159245b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
159347c6ae99SBarry Smith 
159447c6ae99SBarry Smith   PetscFunctionBegin;
159547c6ae99SBarry Smith   /*
159647c6ae99SBarry Smith      nc - number of components per grid point
159747c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
159847c6ae99SBarry Smith   */
15991321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
160047c6ae99SBarry Smith   col  = 2*s + 1;
160147c6ae99SBarry Smith 
1602aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1603aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
160447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
160547c6ae99SBarry Smith 
1606785e854fSJed Brown   ierr = PetscMalloc1(col*col*nc*nc,&cols);CHKERRQ(ierr);
160747c6ae99SBarry Smith 
16081411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
160947c6ae99SBarry Smith 
161047c6ae99SBarry Smith   /* determine the matrix preallocation information */
1611eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
161247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1613bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1614bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
161547c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1616bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1617bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
161847c6ae99SBarry Smith       slot   = i - gxs + gnx*(j - gys);
161947c6ae99SBarry Smith 
162047c6ae99SBarry Smith       /* Find block columns in block row */
162147c6ae99SBarry Smith       cnt = 0;
162247c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
162347c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1624aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) {
162547c6ae99SBarry Smith             cols[cnt++] = slot + ii + gnx*jj;
162647c6ae99SBarry Smith           }
162747c6ae99SBarry Smith         }
162847c6ae99SBarry Smith       }
162945b6f7e9SBarry Smith       ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr);
1630d6e23781SBarry Smith       ierr = MatPreallocateSymmetricSetBlock(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
163147c6ae99SBarry Smith     }
163247c6ae99SBarry Smith   }
163347c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
163447c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
163547c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
163647c6ae99SBarry Smith 
1637784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
163847c6ae99SBarry Smith 
163947c6ae99SBarry Smith   /*
164047c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
164147c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
164247c6ae99SBarry Smith     PETSc ordering.
164347c6ae99SBarry Smith   */
1644fcfd50ebSBarry Smith   if (!da->prealloc_only) {
16451795a4d1SJed Brown     ierr = PetscCalloc1(col*col*nc*nc,&values);CHKERRQ(ierr);
164647c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1647bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1648bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
164947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1650bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1651bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
165247c6ae99SBarry Smith         slot   = i - gxs + gnx*(j - gys);
165347c6ae99SBarry Smith 
165447c6ae99SBarry Smith         /* Find block columns in block row */
165547c6ae99SBarry Smith         cnt = 0;
165647c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
165747c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1658aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) {
165947c6ae99SBarry Smith               cols[cnt++] = slot + ii + gnx*jj;
166047c6ae99SBarry Smith             }
166147c6ae99SBarry Smith           }
166247c6ae99SBarry Smith         }
166345b6f7e9SBarry Smith         ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr);
166447c6ae99SBarry Smith         ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
166547c6ae99SBarry Smith       }
166647c6ae99SBarry Smith     }
166747c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
166847c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166947c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1670189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
167147c6ae99SBarry Smith   }
167247c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
167347c6ae99SBarry Smith   PetscFunctionReturn(0);
167447c6ae99SBarry Smith }
167547c6ae99SBarry Smith 
167647c6ae99SBarry Smith #undef __FUNCT__
1677950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPISBAIJ"
1678950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPISBAIJ(DM da,Mat J)
167947c6ae99SBarry Smith {
168047c6ae99SBarry Smith   PetscErrorCode         ierr;
168147c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
168247c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
168347c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
168447c6ae99SBarry Smith   MPI_Comm               comm;
168547c6ae99SBarry Smith   PetscScalar            *values;
1686bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
1687aa219208SBarry Smith   DMDAStencilType        st;
168845b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
168947c6ae99SBarry Smith 
169047c6ae99SBarry Smith   PetscFunctionBegin;
169147c6ae99SBarry Smith   /*
169247c6ae99SBarry Smith      nc - number of components per grid point
169347c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
169447c6ae99SBarry Smith   */
16951321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
169647c6ae99SBarry Smith   col  = 2*s + 1;
169747c6ae99SBarry Smith 
1698aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1699aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
170047c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
170147c6ae99SBarry Smith 
170247c6ae99SBarry Smith   /* create the matrix */
1703785e854fSJed Brown   ierr = PetscMalloc1(col*col*col,&cols);CHKERRQ(ierr);
170447c6ae99SBarry Smith 
17051411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
170647c6ae99SBarry Smith 
170747c6ae99SBarry Smith   /* determine the matrix preallocation information */
1708eabe889fSLisandro Dalcin   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
170947c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1710bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1711bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
171247c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1713bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1714bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
171547c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1716bff4a2f0SMatthew G. Knepley         kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1717bff4a2f0SMatthew G. Knepley         kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
171847c6ae99SBarry Smith 
171947c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
172047c6ae99SBarry Smith 
172147c6ae99SBarry Smith         /* Find block columns in block row */
172247c6ae99SBarry Smith         cnt = 0;
172347c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
172447c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
172547c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1726aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
172747c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
172847c6ae99SBarry Smith               }
172947c6ae99SBarry Smith             }
173047c6ae99SBarry Smith           }
173147c6ae99SBarry Smith         }
173245b6f7e9SBarry Smith         ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr);
1733d6e23781SBarry Smith         ierr = MatPreallocateSymmetricSetBlock(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
173447c6ae99SBarry Smith       }
173547c6ae99SBarry Smith     }
173647c6ae99SBarry Smith   }
173747c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
173847c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
173947c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
174047c6ae99SBarry Smith 
1741784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
174247c6ae99SBarry Smith 
174347c6ae99SBarry Smith   /*
174447c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
174547c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
174647c6ae99SBarry Smith     PETSc ordering.
174747c6ae99SBarry Smith   */
1748fcfd50ebSBarry Smith   if (!da->prealloc_only) {
17491795a4d1SJed Brown     ierr = PetscCalloc1(col*col*col*nc*nc,&values);CHKERRQ(ierr);
175047c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1751bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1752bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
175347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1754bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1755bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
175647c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
1757bff4a2f0SMatthew G. Knepley           kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1758bff4a2f0SMatthew G. Knepley           kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
175947c6ae99SBarry Smith 
176047c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
176147c6ae99SBarry Smith 
176247c6ae99SBarry Smith           cnt = 0;
176347c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
176447c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
176547c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1766aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
176747c6ae99SBarry Smith                   cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
176847c6ae99SBarry Smith                 }
176947c6ae99SBarry Smith               }
177047c6ae99SBarry Smith             }
177147c6ae99SBarry Smith           }
177245b6f7e9SBarry Smith           ierr = L2GFilterUpperTriangular(ltog,&slot,&cnt,cols);CHKERRQ(ierr);
177347c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
177447c6ae99SBarry Smith         }
177547c6ae99SBarry Smith       }
177647c6ae99SBarry Smith     }
177747c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
177847c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177947c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1780189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
178147c6ae99SBarry Smith   }
178247c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
178347c6ae99SBarry Smith   PetscFunctionReturn(0);
178447c6ae99SBarry Smith }
178547c6ae99SBarry Smith 
178647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
178747c6ae99SBarry Smith 
178847c6ae99SBarry Smith #undef __FUNCT__
1789950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_DA_3d_MPIAIJ_Fill"
1790950540a4SJed Brown PetscErrorCode DMCreateMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J)
179147c6ae99SBarry Smith {
179247c6ae99SBarry Smith   PetscErrorCode         ierr;
179347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1794c0ab637bSBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt, maxcnt = 0,l,p,*dnz,*onz;
179547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
179647c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
179747c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
179847c6ae99SBarry Smith   MPI_Comm               comm;
179947c6ae99SBarry Smith   PetscScalar            *values;
1800bff4a2f0SMatthew G. Knepley   DMBoundaryType         bx,by,bz;
180145b6f7e9SBarry Smith   ISLocalToGlobalMapping ltog;
1802aa219208SBarry Smith   DMDAStencilType        st;
180347c6ae99SBarry Smith 
180447c6ae99SBarry Smith   PetscFunctionBegin;
180547c6ae99SBarry Smith   /*
180647c6ae99SBarry Smith          nc - number of components per grid point
180747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
180847c6ae99SBarry Smith 
180947c6ae99SBarry Smith   */
18101321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
181147c6ae99SBarry Smith   col  = 2*s + 1;
1812bff4a2f0SMatthew G. Knepley   if (bx == DM_BOUNDARY_PERIODIC && (m % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
181347c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1814bff4a2f0SMatthew G. Knepley   if (by == DM_BOUNDARY_PERIODIC && (n % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
181547c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1816bff4a2f0SMatthew G. Knepley   if (bz == DM_BOUNDARY_PERIODIC && (p % col)) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
181747c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
181847c6ae99SBarry Smith 
1819aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1820aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
182147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
182247c6ae99SBarry Smith 
1823785e854fSJed Brown   ierr = PetscMalloc1(col*col*col*nc,&cols);CHKERRQ(ierr);
18241411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
182547c6ae99SBarry Smith 
182647c6ae99SBarry Smith   /* determine the matrix preallocation information */
182747c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
182847c6ae99SBarry Smith 
182906ca8cadSBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
183047c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1831bff4a2f0SMatthew G. Knepley     istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1832bff4a2f0SMatthew G. Knepley     iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
183347c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1834bff4a2f0SMatthew G. Knepley       jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1835bff4a2f0SMatthew G. Knepley       jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
183647c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1837bff4a2f0SMatthew G. Knepley         kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1838bff4a2f0SMatthew G. Knepley         kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
183947c6ae99SBarry Smith 
184047c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
184147c6ae99SBarry Smith 
184247c6ae99SBarry Smith         for (l=0; l<nc; l++) {
184347c6ae99SBarry Smith           cnt = 0;
184447c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
184547c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
184647c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
184747c6ae99SBarry Smith                 if (ii || jj || kk) {
1848aa219208SBarry Smith                   if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
18498865f1eaSKarl 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);
185047c6ae99SBarry Smith                   }
185147c6ae99SBarry Smith                 } else {
185247c6ae99SBarry Smith                   if (dfill) {
18538865f1eaSKarl 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);
185447c6ae99SBarry Smith                   } else {
18558865f1eaSKarl Rupp                     for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
185647c6ae99SBarry Smith                   }
185747c6ae99SBarry Smith                 }
185847c6ae99SBarry Smith               }
185947c6ae99SBarry Smith             }
186047c6ae99SBarry Smith           }
186147c6ae99SBarry Smith           row  = l + nc*(slot);
1862c0ab637bSBarry Smith           maxcnt = PetscMax(maxcnt,cnt);
1863784ac674SJed Brown           ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
186447c6ae99SBarry Smith         }
186547c6ae99SBarry Smith       }
186647c6ae99SBarry Smith     }
186747c6ae99SBarry Smith   }
186847c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
186947c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
187047c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1871784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
187247c6ae99SBarry Smith 
187347c6ae99SBarry Smith   /*
187447c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
187547c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
187647c6ae99SBarry Smith     PETSc ordering.
187747c6ae99SBarry Smith   */
1878fcfd50ebSBarry Smith   if (!da->prealloc_only) {
1879c0ab637bSBarry Smith     ierr = PetscCalloc1(maxcnt,&values);CHKERRQ(ierr);
188047c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1881bff4a2f0SMatthew G. Knepley       istart = (bx == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1882bff4a2f0SMatthew G. Knepley       iend   = (bx == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
188347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1884bff4a2f0SMatthew G. Knepley         jstart = (by == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1885bff4a2f0SMatthew G. Knepley         jend   = (by == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
188647c6ae99SBarry Smith         for (k=zs; k<zs+nz; k++) {
1887bff4a2f0SMatthew G. Knepley           kstart = (bz == DM_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1888bff4a2f0SMatthew G. Knepley           kend   = (bz == DM_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
188947c6ae99SBarry Smith 
189047c6ae99SBarry Smith           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
189147c6ae99SBarry Smith 
189247c6ae99SBarry Smith           for (l=0; l<nc; l++) {
189347c6ae99SBarry Smith             cnt = 0;
189447c6ae99SBarry Smith             for (ii=istart; ii<iend+1; ii++) {
189547c6ae99SBarry Smith               for (jj=jstart; jj<jend+1; jj++) {
189647c6ae99SBarry Smith                 for (kk=kstart; kk<kend+1; kk++) {
189747c6ae99SBarry Smith                   if (ii || jj || kk) {
1898aa219208SBarry Smith                     if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
18998865f1eaSKarl 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);
190047c6ae99SBarry Smith                     }
190147c6ae99SBarry Smith                   } else {
190247c6ae99SBarry Smith                     if (dfill) {
19038865f1eaSKarl 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);
190447c6ae99SBarry Smith                     } else {
19058865f1eaSKarl Rupp                       for (ifill_col=0; ifill_col<nc; ifill_col++) cols[cnt++] = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
190647c6ae99SBarry Smith                     }
190747c6ae99SBarry Smith                   }
190847c6ae99SBarry Smith                 }
190947c6ae99SBarry Smith               }
191047c6ae99SBarry Smith             }
191147c6ae99SBarry Smith             row  = l + nc*(slot);
191247c6ae99SBarry Smith             ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
191347c6ae99SBarry Smith           }
191447c6ae99SBarry Smith         }
191547c6ae99SBarry Smith       }
191647c6ae99SBarry Smith     }
191747c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
191847c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
191947c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1920189e4007SBarry Smith     ierr = MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
192147c6ae99SBarry Smith   }
192247c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
192347c6ae99SBarry Smith   PetscFunctionReturn(0);
192447c6ae99SBarry Smith }
1925