xref: /petsc/src/dm/impls/da/fdda.c (revision 95ee5b0e579985e6d4fa80805eef7b086239d74e)
147c6ae99SBarry Smith 
2c6db04a5SJed Brown #include <private/daimpl.h> /*I      "petscdmda.h"     I*/
3c6db04a5SJed Brown #include <petscmat.h>         /*I      "petscmat.h"    I*/
4c6db04a5SJed Brown #include <private/matimpl.h>
547c6ae99SBarry Smith 
609573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring *);
709573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring *);
809573ac7SBarry Smith extern PetscErrorCode DMGetColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring *);
909573ac7SBarry Smith extern PetscErrorCode DMGetColoring_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"
19aa219208SBarry Smith static PetscErrorCode DMDASetBlockFills_Private(PetscInt *dfill,PetscInt w,PetscInt **rfill)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   PetscErrorCode ierr;
2247c6ae99SBarry Smith   PetscInt       i,j,nz,*fill;
2347c6ae99SBarry Smith 
2447c6ae99SBarry Smith   PetscFunctionBegin;
2547c6ae99SBarry Smith   if (!dfill) PetscFunctionReturn(0);
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   /* count number nonzeros */
2847c6ae99SBarry Smith   nz = 0;
2947c6ae99SBarry Smith   for (i=0; i<w; i++) {
3047c6ae99SBarry Smith     for (j=0; j<w; j++) {
3147c6ae99SBarry Smith       if (dfill[w*i+j]) nz++;
3247c6ae99SBarry Smith     }
3347c6ae99SBarry Smith   }
3447c6ae99SBarry Smith   ierr = PetscMalloc((nz + w + 1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
3547c6ae99SBarry Smith   /* construct modified CSR storage of nonzero structure */
3647c6ae99SBarry Smith   nz = w + 1;
3747c6ae99SBarry Smith   for (i=0; i<w; i++) {
3847c6ae99SBarry Smith     fill[i] = nz;
3947c6ae99SBarry Smith     for (j=0; j<w; j++) {
4047c6ae99SBarry Smith       if (dfill[w*i+j]) {
4147c6ae99SBarry Smith 	fill[nz] = j;
4247c6ae99SBarry Smith 	nz++;
4347c6ae99SBarry Smith       }
4447c6ae99SBarry Smith     }
4547c6ae99SBarry Smith   }
4647c6ae99SBarry Smith   fill[w] = nz;
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith   *rfill = fill;
4947c6ae99SBarry Smith   PetscFunctionReturn(0);
5047c6ae99SBarry Smith }
5147c6ae99SBarry Smith 
5247c6ae99SBarry Smith #undef __FUNCT__
53aa219208SBarry Smith #define __FUNCT__ "DMDASetMatPreallocateOnly"
5447c6ae99SBarry Smith /*@
55aa219208SBarry Smith     DMDASetMatPreallocateOnly - When DMGetMatrix() is called the matrix will be properly
5647c6ae99SBarry Smith        preallocated but the nonzero structure and zero values will not be set.
5747c6ae99SBarry Smith 
58aa219208SBarry Smith     Logically Collective on DMDA
5947c6ae99SBarry Smith 
6047c6ae99SBarry Smith     Input Parameter:
6147c6ae99SBarry Smith +   da - the distributed array
6247c6ae99SBarry Smith -   only - PETSC_TRUE if only want preallocation
6347c6ae99SBarry Smith 
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith     Level: developer
6647c6ae99SBarry Smith 
67aa219208SBarry Smith .seealso DMGetMatrix(), DMDASetGetMatrix(), DMDASetBlockSize(), DMDASetBlockFills()
6847c6ae99SBarry Smith 
6947c6ae99SBarry Smith @*/
707087cfbeSBarry Smith PetscErrorCode  DMDASetMatPreallocateOnly(DM da,PetscBool  only)
7147c6ae99SBarry Smith {
7247c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
7347c6ae99SBarry Smith 
7447c6ae99SBarry Smith   PetscFunctionBegin;
7547c6ae99SBarry Smith   dd->prealloc_only = only;
7647c6ae99SBarry Smith   PetscFunctionReturn(0);
7747c6ae99SBarry Smith }
7847c6ae99SBarry Smith 
7947c6ae99SBarry Smith #undef __FUNCT__
80aa219208SBarry Smith #define __FUNCT__ "DMDASetBlockFills"
8147c6ae99SBarry Smith /*@
82aa219208SBarry Smith     DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem
8394013140SBarry Smith     of the matrix returned by DMGetMatrix().
8447c6ae99SBarry Smith 
85aa219208SBarry Smith     Logically Collective on DMDA
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith     Input Parameter:
8847c6ae99SBarry Smith +   da - the distributed array
8947c6ae99SBarry Smith .   dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block)
9047c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
9147c6ae99SBarry Smith 
9247c6ae99SBarry Smith 
9347c6ae99SBarry Smith     Level: developer
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith     Notes: This only makes sense when you are doing multicomponent problems but using the
9647c6ae99SBarry Smith        MPIAIJ matrix format
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
9947c6ae99SBarry Smith        representing coupling and 0 entries for missing coupling. For example
10047c6ae99SBarry Smith $             dfill[9] = {1, 0, 0,
10147c6ae99SBarry Smith $                         1, 1, 0,
10247c6ae99SBarry Smith $                         0, 1, 1}
10347c6ae99SBarry Smith        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
10447c6ae99SBarry Smith        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
10547c6ae99SBarry Smith        diagonal block).
10647c6ae99SBarry Smith 
107aa219208SBarry Smith      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
10847c6ae99SBarry Smith      can be represented in the dfill, ofill format
10947c6ae99SBarry Smith 
11047c6ae99SBarry Smith    Contributed by Glenn Hammond
11147c6ae99SBarry Smith 
112aa219208SBarry Smith .seealso DMGetMatrix(), DMDASetGetMatrix(), DMDASetMatPreallocateOnly()
11347c6ae99SBarry Smith 
11447c6ae99SBarry Smith @*/
1157087cfbeSBarry Smith PetscErrorCode  DMDASetBlockFills(DM da,PetscInt *dfill,PetscInt *ofill)
11647c6ae99SBarry Smith {
11747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
11847c6ae99SBarry Smith   PetscErrorCode ierr;
11947c6ae99SBarry Smith 
12047c6ae99SBarry Smith   PetscFunctionBegin;
121aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr);
122aa219208SBarry Smith   ierr = DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr);
12347c6ae99SBarry Smith   PetscFunctionReturn(0);
12447c6ae99SBarry Smith }
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith #undef __FUNCT__
12894013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA"
1297087cfbeSBarry Smith PetscErrorCode  DMGetColoring_DA(DM da,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
13047c6ae99SBarry Smith {
13147c6ae99SBarry Smith   PetscErrorCode   ierr;
13247c6ae99SBarry Smith   PetscInt         dim,m,n,p,nc;
1331321219cSEthan Coon   DMDABoundaryType bx,by,bz;
13447c6ae99SBarry Smith   MPI_Comm         comm;
13547c6ae99SBarry Smith   PetscMPIInt      size;
13647c6ae99SBarry Smith   PetscBool        isBAIJ;
13747c6ae99SBarry Smith   DM_DA            *dd = (DM_DA*)da->data;
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith   PetscFunctionBegin;
14047c6ae99SBarry Smith   /*
14147c6ae99SBarry Smith                                   m
14247c6ae99SBarry Smith           ------------------------------------------------------
14347c6ae99SBarry Smith          |                                                     |
14447c6ae99SBarry Smith          |                                                     |
14547c6ae99SBarry Smith          |               ----------------------                |
14647c6ae99SBarry Smith          |               |                    |                |
14747c6ae99SBarry Smith       n  |           yn  |                    |                |
14847c6ae99SBarry Smith          |               |                    |                |
14947c6ae99SBarry Smith          |               .---------------------                |
15047c6ae99SBarry Smith          |             (xs,ys)     xn                          |
15147c6ae99SBarry Smith          |            .                                        |
15247c6ae99SBarry Smith          |         (gxs,gys)                                   |
15347c6ae99SBarry Smith          |                                                     |
15447c6ae99SBarry Smith           -----------------------------------------------------
15547c6ae99SBarry Smith   */
15647c6ae99SBarry Smith 
15747c6ae99SBarry Smith   /*
15847c6ae99SBarry Smith          nc - number of components per grid point
15947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
16047c6ae99SBarry Smith 
16147c6ae99SBarry Smith   */
1621321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&bx,&by,&bz,0);CHKERRQ(ierr);
16347c6ae99SBarry Smith 
16447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
16547c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
16647c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED){
16747c6ae99SBarry Smith     if (size == 1) {
16847c6ae99SBarry Smith       ctype = IS_COLORING_GLOBAL;
16947c6ae99SBarry Smith     } else if (dim > 1){
1701321219cSEthan Coon       if ((m==1 && bx == DMDA_BOUNDARY_PERIODIC) || (n==1 && by == DMDA_BOUNDARY_PERIODIC) || (p==1 && bz == DMDA_BOUNDARY_PERIODIC)){
17147c6ae99SBarry Smith         SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"IS_COLORING_GHOSTED cannot be used for periodic boundary condition having both ends of the domain  on the same process");
17247c6ae99SBarry Smith       }
17347c6ae99SBarry Smith     }
17447c6ae99SBarry Smith   }
17547c6ae99SBarry Smith 
176aa219208SBarry Smith   /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ
17747c6ae99SBarry Smith      matrices is for the blocks, not the individual matrix elements  */
17847c6ae99SBarry Smith   ierr = PetscStrcmp(mtype,MATBAIJ,&isBAIJ);CHKERRQ(ierr);
17947c6ae99SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);}
18047c6ae99SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);}
18147c6ae99SBarry Smith   if (isBAIJ) {
18247c6ae99SBarry Smith     dd->w = 1;
18347c6ae99SBarry Smith     dd->xs = dd->xs/nc;
18447c6ae99SBarry Smith     dd->xe = dd->xe/nc;
18547c6ae99SBarry Smith     dd->Xs = dd->Xs/nc;
18647c6ae99SBarry Smith     dd->Xe = dd->Xe/nc;
18747c6ae99SBarry Smith   }
18847c6ae99SBarry Smith 
18947c6ae99SBarry Smith   /*
190aa219208SBarry Smith      We do not provide a getcoloring function in the DMDA operations because
191aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
19247c6ae99SBarry Smith    more low-level then matrices.
19347c6ae99SBarry Smith   */
19447c6ae99SBarry Smith   if (dim == 1) {
19594013140SBarry Smith     ierr = DMGetColoring_DA_1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
19647c6ae99SBarry Smith   } else if (dim == 2) {
19794013140SBarry Smith     ierr =  DMGetColoring_DA_2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
19847c6ae99SBarry Smith   } else if (dim == 3) {
19994013140SBarry Smith     ierr =  DMGetColoring_DA_3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
20047c6ae99SBarry Smith   } else {
20147c6ae99SBarry Smith     SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim);
20247c6ae99SBarry Smith   }
20347c6ae99SBarry Smith   if (isBAIJ) {
20447c6ae99SBarry Smith     dd->w = nc;
20547c6ae99SBarry Smith     dd->xs = dd->xs*nc;
20647c6ae99SBarry Smith     dd->xe = dd->xe*nc;
20747c6ae99SBarry Smith     dd->Xs = dd->Xs*nc;
20847c6ae99SBarry Smith     dd->Xe = dd->Xe*nc;
20947c6ae99SBarry Smith   }
21047c6ae99SBarry Smith   PetscFunctionReturn(0);
21147c6ae99SBarry Smith }
21247c6ae99SBarry Smith 
21347c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith #undef __FUNCT__
21694013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_2d_MPIAIJ"
21794013140SBarry Smith PetscErrorCode DMGetColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
21847c6ae99SBarry Smith {
21947c6ae99SBarry Smith   PetscErrorCode         ierr;
22047c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col;
22147c6ae99SBarry Smith   PetscInt               ncolors;
22247c6ae99SBarry Smith   MPI_Comm               comm;
2231321219cSEthan Coon   DMDABoundaryType       bx,by;
224aa219208SBarry Smith   DMDAStencilType        st;
22547c6ae99SBarry Smith   ISColoringValue        *colors;
22647c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
22747c6ae99SBarry Smith 
22847c6ae99SBarry Smith   PetscFunctionBegin;
22947c6ae99SBarry Smith   /*
23047c6ae99SBarry Smith          nc - number of components per grid point
23147c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
23247c6ae99SBarry Smith 
23347c6ae99SBarry Smith   */
2341321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
23547c6ae99SBarry Smith   col    = 2*s + 1;
236aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
237aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
23847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
23947c6ae99SBarry Smith 
24047c6ae99SBarry Smith   /* special case as taught to us by Paul Hovland */
241aa219208SBarry Smith   if (st == DMDA_STENCIL_STAR && s == 1) {
24294013140SBarry Smith     ierr = DMGetColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
24347c6ae99SBarry Smith   } else {
24447c6ae99SBarry Smith 
2451321219cSEthan Coon     if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){
24647c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\
24747c6ae99SBarry Smith                  by 2*stencil_width + 1 (%d)\n", m, col);
24847c6ae99SBarry Smith     }
2491321219cSEthan Coon     if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){
25047c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\
25147c6ae99SBarry Smith                  by 2*stencil_width + 1 (%d)\n", n, col);
25247c6ae99SBarry Smith     }
25347c6ae99SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
25447c6ae99SBarry Smith       if (!dd->localcoloring) {
25547c6ae99SBarry Smith 	ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
25647c6ae99SBarry Smith 	ii = 0;
25747c6ae99SBarry Smith 	for (j=ys; j<ys+ny; j++) {
25847c6ae99SBarry Smith 	  for (i=xs; i<xs+nx; i++) {
25947c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
26047c6ae99SBarry Smith 	      colors[ii++] = k + nc*((i % col) + col*(j % col));
26147c6ae99SBarry Smith 	    }
26247c6ae99SBarry Smith 	  }
26347c6ae99SBarry Smith 	}
26447c6ae99SBarry Smith         ncolors = nc + nc*(col-1 + col*(col-1));
26547c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
26647c6ae99SBarry Smith       }
26747c6ae99SBarry Smith       *coloring = dd->localcoloring;
26847c6ae99SBarry Smith     } else if (ctype == IS_COLORING_GHOSTED) {
26947c6ae99SBarry Smith       if (!dd->ghostedcoloring) {
27047c6ae99SBarry Smith 	ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
27147c6ae99SBarry Smith 	ii = 0;
27247c6ae99SBarry Smith 	for (j=gys; j<gys+gny; j++) {
27347c6ae99SBarry Smith 	  for (i=gxs; i<gxs+gnx; i++) {
27447c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
27547c6ae99SBarry Smith 	      /* the complicated stuff is to handle periodic boundaries */
27647c6ae99SBarry Smith 	      colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
27747c6ae99SBarry Smith 	    }
27847c6ae99SBarry Smith 	  }
27947c6ae99SBarry Smith 	}
28047c6ae99SBarry Smith         ncolors = nc + nc*(col - 1 + col*(col-1));
28147c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
28247c6ae99SBarry Smith         /* PetscIntView(ncolors,(PetscInt *)colors,0); */
28347c6ae99SBarry Smith 
28447c6ae99SBarry Smith 	ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
28547c6ae99SBarry Smith       }
28647c6ae99SBarry Smith       *coloring = dd->ghostedcoloring;
28747c6ae99SBarry Smith     } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
28847c6ae99SBarry Smith   }
28947c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
29047c6ae99SBarry Smith   PetscFunctionReturn(0);
29147c6ae99SBarry Smith }
29247c6ae99SBarry Smith 
29347c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
29447c6ae99SBarry Smith 
29547c6ae99SBarry Smith #undef __FUNCT__
29694013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_3d_MPIAIJ"
29794013140SBarry Smith PetscErrorCode DMGetColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
29847c6ae99SBarry Smith {
29947c6ae99SBarry Smith   PetscErrorCode    ierr;
30047c6ae99SBarry 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;
30147c6ae99SBarry Smith   PetscInt          ncolors;
30247c6ae99SBarry Smith   MPI_Comm          comm;
3031321219cSEthan Coon   DMDABoundaryType  bx,by,bz;
304aa219208SBarry Smith   DMDAStencilType   st;
30547c6ae99SBarry Smith   ISColoringValue   *colors;
30647c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
30747c6ae99SBarry Smith 
30847c6ae99SBarry Smith   PetscFunctionBegin;
30947c6ae99SBarry Smith   /*
31047c6ae99SBarry Smith          nc - number of components per grid point
31147c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
31247c6ae99SBarry Smith 
31347c6ae99SBarry Smith   */
3141321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
31547c6ae99SBarry Smith   col    = 2*s + 1;
3161321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){
31747c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
31847c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
31947c6ae99SBarry Smith   }
3201321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){
32147c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
32247c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
32347c6ae99SBarry Smith   }
3241321219cSEthan Coon   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)){
32547c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
32647c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
32747c6ae99SBarry Smith   }
32847c6ae99SBarry Smith 
329aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
330aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
33147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
33247c6ae99SBarry Smith 
33347c6ae99SBarry Smith   /* create the coloring */
33447c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
33547c6ae99SBarry Smith     if (!dd->localcoloring) {
33647c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
33747c6ae99SBarry Smith       ii = 0;
33847c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
33947c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
34047c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
34147c6ae99SBarry Smith             for (l=0; l<nc; l++) {
34247c6ae99SBarry Smith               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
34347c6ae99SBarry Smith             }
34447c6ae99SBarry Smith           }
34547c6ae99SBarry Smith         }
34647c6ae99SBarry Smith       }
34747c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
34847c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr);
34947c6ae99SBarry Smith     }
35047c6ae99SBarry Smith     *coloring = dd->localcoloring;
35147c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
35247c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
35347c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
35447c6ae99SBarry Smith       ii = 0;
35547c6ae99SBarry Smith       for (k=gzs; k<gzs+gnz; k++) {
35647c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
35747c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
35847c6ae99SBarry Smith             for (l=0; l<nc; l++) {
35947c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
36047c6ae99SBarry Smith               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
36147c6ae99SBarry Smith             }
36247c6ae99SBarry Smith           }
36347c6ae99SBarry Smith         }
36447c6ae99SBarry Smith       }
36547c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
36647c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
36747c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
36847c6ae99SBarry Smith     }
36947c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
37047c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
37147c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
37247c6ae99SBarry Smith   PetscFunctionReturn(0);
37347c6ae99SBarry Smith }
37447c6ae99SBarry Smith 
37547c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
37647c6ae99SBarry Smith 
37747c6ae99SBarry Smith #undef __FUNCT__
37894013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_1d_MPIAIJ"
37994013140SBarry Smith PetscErrorCode DMGetColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
38047c6ae99SBarry Smith {
38147c6ae99SBarry Smith   PetscErrorCode    ierr;
38247c6ae99SBarry Smith   PetscInt          xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
38347c6ae99SBarry Smith   PetscInt          ncolors;
38447c6ae99SBarry Smith   MPI_Comm          comm;
3851321219cSEthan Coon   DMDABoundaryType  bx;
38647c6ae99SBarry Smith   ISColoringValue   *colors;
38747c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
38847c6ae99SBarry Smith 
38947c6ae99SBarry Smith   PetscFunctionBegin;
39047c6ae99SBarry Smith   /*
39147c6ae99SBarry Smith          nc - number of components per grid point
39247c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
39347c6ae99SBarry Smith 
39447c6ae99SBarry Smith   */
3951321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
39647c6ae99SBarry Smith   col    = 2*s + 1;
39747c6ae99SBarry Smith 
3981321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) {
39947c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points is divisible\n\
40047c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
40147c6ae99SBarry Smith   }
40247c6ae99SBarry Smith 
403aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
404aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
40547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
40647c6ae99SBarry Smith 
40747c6ae99SBarry Smith   /* create the coloring */
40847c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
40947c6ae99SBarry Smith     if (!dd->localcoloring) {
41047c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
41147c6ae99SBarry Smith       i1 = 0;
41247c6ae99SBarry Smith       for (i=xs; i<xs+nx; i++) {
41347c6ae99SBarry Smith         for (l=0; l<nc; l++) {
41447c6ae99SBarry Smith           colors[i1++] = l + nc*(i % col);
41547c6ae99SBarry Smith         }
41647c6ae99SBarry Smith       }
41747c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
41847c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr);
41947c6ae99SBarry Smith     }
42047c6ae99SBarry Smith     *coloring = dd->localcoloring;
42147c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
42247c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
42347c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
42447c6ae99SBarry Smith       i1 = 0;
42547c6ae99SBarry Smith       for (i=gxs; i<gxs+gnx; i++) {
42647c6ae99SBarry Smith         for (l=0; l<nc; l++) {
42747c6ae99SBarry Smith           /* the complicated stuff is to handle periodic boundaries */
42847c6ae99SBarry Smith           colors[i1++] = l + nc*(SetInRange(i,m) % col);
42947c6ae99SBarry Smith         }
43047c6ae99SBarry Smith       }
43147c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
43247c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
43347c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
43447c6ae99SBarry Smith     }
43547c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
43647c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
43747c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
43847c6ae99SBarry Smith   PetscFunctionReturn(0);
43947c6ae99SBarry Smith }
44047c6ae99SBarry Smith 
44147c6ae99SBarry Smith #undef __FUNCT__
44294013140SBarry Smith #define __FUNCT__ "DMGetColoring_DA_2d_5pt_MPIAIJ"
44394013140SBarry Smith PetscErrorCode DMGetColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
44447c6ae99SBarry Smith {
44547c6ae99SBarry Smith   PetscErrorCode    ierr;
44647c6ae99SBarry Smith   PetscInt          xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
44747c6ae99SBarry Smith   PetscInt          ncolors;
44847c6ae99SBarry Smith   MPI_Comm          comm;
4491321219cSEthan Coon   DMDABoundaryType  bx,by;
45047c6ae99SBarry Smith   ISColoringValue   *colors;
45147c6ae99SBarry Smith   DM_DA             *dd = (DM_DA*)da->data;
45247c6ae99SBarry Smith 
45347c6ae99SBarry Smith   PetscFunctionBegin;
45447c6ae99SBarry Smith   /*
45547c6ae99SBarry Smith          nc - number of components per grid point
45647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
45747c6ae99SBarry Smith 
45847c6ae99SBarry Smith   */
4591321219cSEthan Coon   ierr   = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);CHKERRQ(ierr);
460aa219208SBarry Smith   ierr   = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
461aa219208SBarry Smith   ierr   = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
46247c6ae99SBarry Smith   ierr   = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
46347c6ae99SBarry Smith 
4641321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && (m % 5)){
46547c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
46647c6ae99SBarry Smith                  by 5\n");
46747c6ae99SBarry Smith   }
4681321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC && (n % 5)){
46947c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
47047c6ae99SBarry Smith                  by 5\n");
47147c6ae99SBarry Smith   }
47247c6ae99SBarry Smith 
47347c6ae99SBarry Smith   /* create the coloring */
47447c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
47547c6ae99SBarry Smith     if (!dd->localcoloring) {
47647c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
47747c6ae99SBarry Smith       ii = 0;
47847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
47947c6ae99SBarry Smith 	for (i=xs; i<xs+nx; i++) {
48047c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
48147c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*j+i) % 5);
48247c6ae99SBarry Smith 	  }
48347c6ae99SBarry Smith 	}
48447c6ae99SBarry Smith       }
48547c6ae99SBarry Smith       ncolors = 5*nc;
48647c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
48747c6ae99SBarry Smith     }
48847c6ae99SBarry Smith     *coloring = dd->localcoloring;
48947c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
49047c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
49147c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
49247c6ae99SBarry Smith       ii = 0;
49347c6ae99SBarry Smith       for (j=gys; j<gys+gny; j++) {
49447c6ae99SBarry Smith 	for (i=gxs; i<gxs+gnx; i++) {
49547c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
49647c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
49747c6ae99SBarry Smith 	  }
49847c6ae99SBarry Smith 	}
49947c6ae99SBarry Smith       }
50047c6ae99SBarry Smith       ncolors = 5*nc;
50147c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
50247c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
50347c6ae99SBarry Smith     }
50447c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
50547c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
50647c6ae99SBarry Smith   PetscFunctionReturn(0);
50747c6ae99SBarry Smith }
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith /* =========================================================================== */
51009573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_1d_MPIAIJ(DM,Mat);
51109573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ(DM,Mat);
51209573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ_Fill(DM,Mat);
51309573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ(DM,Mat);
51409573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ_Fill(DM,Mat);
51509573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPIBAIJ(DM,Mat);
51609573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPIBAIJ(DM,Mat);
51709573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_2d_MPISBAIJ(DM,Mat);
51809573ac7SBarry Smith extern PetscErrorCode DMGetMatrix_DA_3d_MPISBAIJ(DM,Mat);
51947c6ae99SBarry Smith 
52047c6ae99SBarry Smith #undef __FUNCT__
521*95ee5b0eSBarry Smith #define __FUNCT__ "MatSetDM"
52247c6ae99SBarry Smith /*@
523*95ee5b0eSBarry Smith    MatSetDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix
52447c6ae99SBarry Smith 
52547c6ae99SBarry Smith    Logically Collective on Mat
52647c6ae99SBarry Smith 
52747c6ae99SBarry Smith    Input Parameters:
52847c6ae99SBarry Smith +  mat - the matrix
52947c6ae99SBarry Smith -  da - the da
53047c6ae99SBarry Smith 
53147c6ae99SBarry Smith    Level: intermediate
53247c6ae99SBarry Smith 
53347c6ae99SBarry Smith @*/
534*95ee5b0eSBarry Smith PetscErrorCode  MatSetDM(Mat mat,DM da)
53547c6ae99SBarry Smith {
53647c6ae99SBarry Smith   PetscErrorCode ierr;
53747c6ae99SBarry Smith 
53847c6ae99SBarry Smith   PetscFunctionBegin;
53947c6ae99SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
54047c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
541*95ee5b0eSBarry Smith   ierr = PetscTryMethod(mat,"MatSetDM_C",(Mat,DM),(mat,da));CHKERRQ(ierr);
54247c6ae99SBarry Smith   PetscFunctionReturn(0);
54347c6ae99SBarry Smith }
54447c6ae99SBarry Smith 
54547c6ae99SBarry Smith EXTERN_C_BEGIN
54647c6ae99SBarry Smith #undef __FUNCT__
54747c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA"
5487087cfbeSBarry Smith PetscErrorCode  MatView_MPI_DA(Mat A,PetscViewer viewer)
54947c6ae99SBarry Smith {
5509a42bb27SBarry Smith   DM             da;
55147c6ae99SBarry Smith   PetscErrorCode ierr;
55247c6ae99SBarry Smith   const char     *prefix;
55347c6ae99SBarry Smith   Mat            Anatural;
55447c6ae99SBarry Smith   AO             ao;
55547c6ae99SBarry Smith   PetscInt       rstart,rend,*petsc,i;
55647c6ae99SBarry Smith   IS             is;
55747c6ae99SBarry Smith   MPI_Comm       comm;
55847c6ae99SBarry Smith 
55947c6ae99SBarry Smith   PetscFunctionBegin;
56047c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
561aa219208SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"DMDA",(PetscObject*)&da);CHKERRQ(ierr);
562aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
56347c6ae99SBarry Smith 
564aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
56547c6ae99SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
56647c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr);
56747c6ae99SBarry Smith   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
56847c6ae99SBarry Smith   ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr);
56947c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
57047c6ae99SBarry Smith 
57147c6ae99SBarry Smith   /* call viewer on natural ordering */
57247c6ae99SBarry Smith   ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr);
57347c6ae99SBarry Smith   ierr = ISDestroy(is);CHKERRQ(ierr);
57447c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr);
57547c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr);
57647c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr);
57747c6ae99SBarry Smith   ierr = MatView(Anatural,viewer);CHKERRQ(ierr);
57847c6ae99SBarry Smith   ierr = MatDestroy(Anatural);CHKERRQ(ierr);
57947c6ae99SBarry Smith   PetscFunctionReturn(0);
58047c6ae99SBarry Smith }
58147c6ae99SBarry Smith EXTERN_C_END
58247c6ae99SBarry Smith 
58347c6ae99SBarry Smith EXTERN_C_BEGIN
58447c6ae99SBarry Smith #undef __FUNCT__
58547c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA"
5867087cfbeSBarry Smith PetscErrorCode  MatLoad_MPI_DA(Mat A,PetscViewer viewer)
58747c6ae99SBarry Smith {
5889a42bb27SBarry Smith   DM             da;
58947c6ae99SBarry Smith   PetscErrorCode ierr;
59047c6ae99SBarry Smith   Mat            Anatural,Aapp;
59147c6ae99SBarry Smith   AO             ao;
59247c6ae99SBarry Smith   PetscInt       rstart,rend,*app,i;
59347c6ae99SBarry Smith   IS             is;
59447c6ae99SBarry Smith   MPI_Comm       comm;
59547c6ae99SBarry Smith 
59647c6ae99SBarry Smith   PetscFunctionBegin;
59747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
598aa219208SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"DMDA",(PetscObject*)&da);CHKERRQ(ierr);
599aa219208SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");
60047c6ae99SBarry Smith 
60147c6ae99SBarry Smith   /* Load the matrix in natural ordering */
60247c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,&Anatural);CHKERRQ(ierr);
60347c6ae99SBarry Smith   ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr);
60447c6ae99SBarry Smith   ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
60547c6ae99SBarry Smith   ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr);
60647c6ae99SBarry Smith 
60747c6ae99SBarry Smith   /* Map natural ordering to application ordering and create IS */
608aa219208SBarry Smith   ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
60947c6ae99SBarry Smith   ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr);
61047c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr);
61147c6ae99SBarry Smith   for (i=rstart; i<rend; i++) app[i-rstart] = i;
61247c6ae99SBarry Smith   ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr);
61347c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
61447c6ae99SBarry Smith 
61547c6ae99SBarry Smith   /* Do permutation and replace header */
61647c6ae99SBarry Smith   ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr);
61747c6ae99SBarry Smith   ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr);
61847c6ae99SBarry Smith   ierr = ISDestroy(is);CHKERRQ(ierr);
61947c6ae99SBarry Smith   ierr = MatDestroy(Anatural);CHKERRQ(ierr);
62047c6ae99SBarry Smith   PetscFunctionReturn(0);
62147c6ae99SBarry Smith }
62247c6ae99SBarry Smith EXTERN_C_END
62347c6ae99SBarry Smith 
62447c6ae99SBarry Smith #undef __FUNCT__
62594013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA"
6267087cfbeSBarry Smith PetscErrorCode  DMGetMatrix_DA(DM da, const MatType mtype,Mat *J)
62747c6ae99SBarry Smith {
62847c6ae99SBarry Smith   PetscErrorCode ierr;
62947c6ae99SBarry Smith   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
63047c6ae99SBarry Smith   Mat            A;
63147c6ae99SBarry Smith   MPI_Comm       comm;
63247c6ae99SBarry Smith   const MatType  Atype;
63347c6ae99SBarry Smith   void           (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL;
63447c6ae99SBarry Smith   MatType        ttype[256];
63547c6ae99SBarry Smith   PetscBool      flg;
63647c6ae99SBarry Smith   PetscMPIInt    size;
63747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
63847c6ae99SBarry Smith 
63947c6ae99SBarry Smith   PetscFunctionBegin;
64047c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
64147c6ae99SBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
64247c6ae99SBarry Smith #endif
6435da5aae0SJed Brown   if (!mtype) mtype = MATAIJ;
64447c6ae99SBarry Smith   ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr);
645aa219208SBarry Smith   ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DMDA options","Mat");CHKERRQ(ierr);
64647c6ae99SBarry Smith   ierr = PetscOptionsList("-da_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr);
64747c6ae99SBarry Smith   ierr = PetscOptionsEnd();
64847c6ae99SBarry Smith 
64947c6ae99SBarry Smith   /*
65047c6ae99SBarry Smith                                   m
65147c6ae99SBarry Smith           ------------------------------------------------------
65247c6ae99SBarry Smith          |                                                     |
65347c6ae99SBarry Smith          |                                                     |
65447c6ae99SBarry Smith          |               ----------------------                |
65547c6ae99SBarry Smith          |               |                    |                |
65647c6ae99SBarry Smith       n  |           ny  |                    |                |
65747c6ae99SBarry Smith          |               |                    |                |
65847c6ae99SBarry Smith          |               .---------------------                |
65947c6ae99SBarry Smith          |             (xs,ys)     nx                          |
66047c6ae99SBarry Smith          |            .                                        |
66147c6ae99SBarry Smith          |         (gxs,gys)                                   |
66247c6ae99SBarry Smith          |                                                     |
66347c6ae99SBarry Smith           -----------------------------------------------------
66447c6ae99SBarry Smith   */
66547c6ae99SBarry Smith 
66647c6ae99SBarry Smith   /*
66747c6ae99SBarry Smith          nc - number of components per grid point
66847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
66947c6ae99SBarry Smith 
67047c6ae99SBarry Smith   */
6711321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
672aa219208SBarry Smith   ierr = DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
67347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
67447c6ae99SBarry Smith   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
67547c6ae99SBarry Smith   ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr);
67647c6ae99SBarry Smith   ierr = MatSetType(A,(const MatType)ttype);CHKERRQ(ierr);
677*95ee5b0eSBarry Smith   ierr = MatSetDM(A,da);CHKERRQ(ierr);
67847c6ae99SBarry Smith   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
67947c6ae99SBarry Smith   ierr = MatGetType(A,&Atype);CHKERRQ(ierr);
68047c6ae99SBarry Smith   /*
681aa219208SBarry Smith      We do not provide a getmatrix function in the DMDA operations because
682aa219208SBarry Smith    the basic DMDA does not know about matrices. We think of DMDA as being more
68347c6ae99SBarry Smith    more low-level than matrices. This is kind of cheating but, cause sometimes
684aa219208SBarry Smith    we think of DMDA has higher level than matrices.
68547c6ae99SBarry Smith 
68647c6ae99SBarry Smith      We could switch based on Atype (or mtype), but we do not since the
68747c6ae99SBarry Smith    specialized setting routines depend only the particular preallocation
68847c6ae99SBarry Smith    details of the matrix, not the type itself.
68947c6ae99SBarry Smith   */
69047c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
69147c6ae99SBarry Smith   if (!aij) {
69247c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
69347c6ae99SBarry Smith   }
69447c6ae99SBarry Smith   if (!aij) {
69547c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
69647c6ae99SBarry Smith     if (!baij) {
69747c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
69847c6ae99SBarry Smith     }
69947c6ae99SBarry Smith     if (!baij){
70047c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
70147c6ae99SBarry Smith       if (!sbaij) {
70247c6ae99SBarry Smith         ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
70347c6ae99SBarry Smith       }
70447c6ae99SBarry Smith     }
70547c6ae99SBarry Smith   }
70647c6ae99SBarry Smith   if (aij) {
70747c6ae99SBarry Smith     if (dim == 1) {
70894013140SBarry Smith       ierr = DMGetMatrix_DA_1d_MPIAIJ(da,A);CHKERRQ(ierr);
70947c6ae99SBarry Smith     } else if (dim == 2) {
71047c6ae99SBarry Smith       if (dd->ofill) {
71194013140SBarry Smith         ierr = DMGetMatrix_DA_2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
71247c6ae99SBarry Smith       } else {
71394013140SBarry Smith         ierr = DMGetMatrix_DA_2d_MPIAIJ(da,A);CHKERRQ(ierr);
71447c6ae99SBarry Smith       }
71547c6ae99SBarry Smith     } else if (dim == 3) {
71647c6ae99SBarry Smith       if (dd->ofill) {
71794013140SBarry Smith         ierr = DMGetMatrix_DA_3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
71847c6ae99SBarry Smith       } else {
71994013140SBarry Smith         ierr = DMGetMatrix_DA_3d_MPIAIJ(da,A);CHKERRQ(ierr);
72047c6ae99SBarry Smith       }
72147c6ae99SBarry Smith     }
72247c6ae99SBarry Smith   } else if (baij) {
72347c6ae99SBarry Smith     if (dim == 2) {
72494013140SBarry Smith       ierr = DMGetMatrix_DA_2d_MPIBAIJ(da,A);CHKERRQ(ierr);
72547c6ae99SBarry Smith     } else if (dim == 3) {
72694013140SBarry Smith       ierr = DMGetMatrix_DA_3d_MPIBAIJ(da,A);CHKERRQ(ierr);
72747c6ae99SBarry Smith     } else {
728b17742caSSean Farley       SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \
729b17742caSSean Farley 	       "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim);
73047c6ae99SBarry Smith     }
73147c6ae99SBarry Smith   } else if (sbaij) {
73247c6ae99SBarry Smith     if (dim == 2) {
73394013140SBarry Smith       ierr = DMGetMatrix_DA_2d_MPISBAIJ(da,A);CHKERRQ(ierr);
73447c6ae99SBarry Smith     } else if (dim == 3) {
73594013140SBarry Smith       ierr = DMGetMatrix_DA_3d_MPISBAIJ(da,A);CHKERRQ(ierr);
73647c6ae99SBarry Smith     } else {
737b17742caSSean Farley       SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \
738b17742caSSean Farley 	       "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim);
73947c6ae99SBarry Smith     }
74047c6ae99SBarry Smith   }
741aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
74247c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
743aa219208SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"DMDA",(PetscObject)da);CHKERRQ(ierr);
74447c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
74547c6ae99SBarry Smith   if (size > 1) {
74647c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
74747c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);CHKERRQ(ierr);
74847c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);CHKERRQ(ierr);
74947c6ae99SBarry Smith   }
75047c6ae99SBarry Smith   *J = A;
75147c6ae99SBarry Smith   PetscFunctionReturn(0);
75247c6ae99SBarry Smith }
75347c6ae99SBarry Smith 
75447c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
75547c6ae99SBarry Smith #undef __FUNCT__
75694013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIAIJ"
75794013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ(DM da,Mat J)
75847c6ae99SBarry Smith {
75947c6ae99SBarry Smith   PetscErrorCode         ierr;
76047c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny,m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p;
76147c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
76247c6ae99SBarry Smith   MPI_Comm               comm;
76347c6ae99SBarry Smith   PetscScalar            *values;
7641321219cSEthan Coon   DMDABoundaryType       bx,by;
76547c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
766aa219208SBarry Smith   DMDAStencilType        st;
76747c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
76847c6ae99SBarry Smith 
76947c6ae99SBarry Smith   PetscFunctionBegin;
77047c6ae99SBarry Smith   /*
77147c6ae99SBarry Smith          nc - number of components per grid point
77247c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
77347c6ae99SBarry Smith 
77447c6ae99SBarry Smith   */
7751321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
77647c6ae99SBarry Smith   col = 2*s + 1;
777aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
778aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
77947c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
78047c6ae99SBarry Smith 
78147c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
7821411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
7831411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
78447c6ae99SBarry Smith 
78547c6ae99SBarry Smith   /* determine the matrix preallocation information */
78647c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
78747c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
78847c6ae99SBarry Smith 
7891321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
7901321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
79147c6ae99SBarry Smith 
79247c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
79347c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
79447c6ae99SBarry Smith 
7951321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
7961321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
79747c6ae99SBarry Smith 
79847c6ae99SBarry Smith       cnt  = 0;
79947c6ae99SBarry Smith       for (k=0; k<nc; k++) {
80047c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
80147c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
802aa219208SBarry Smith 	    if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
80347c6ae99SBarry Smith 	      cols[cnt++]  = k + nc*(slot + gnx*l + p);
80447c6ae99SBarry Smith 	    }
80547c6ae99SBarry Smith 	  }
80647c6ae99SBarry Smith 	}
80747c6ae99SBarry Smith 	rows[k] = k + nc*(slot);
80847c6ae99SBarry Smith       }
809784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
81047c6ae99SBarry Smith     }
81147c6ae99SBarry Smith   }
81247c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
81347c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
81447c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
81547c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
81647c6ae99SBarry Smith 
817784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
818784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
81947c6ae99SBarry Smith 
82047c6ae99SBarry Smith   /*
82147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
82247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
82347c6ae99SBarry Smith     PETSc ordering.
82447c6ae99SBarry Smith   */
82547c6ae99SBarry Smith   if (!dd->prealloc_only) {
82647c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
82747c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
82847c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
82947c6ae99SBarry Smith 
8301321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8311321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
83247c6ae99SBarry Smith 
83347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
83447c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
83547c6ae99SBarry Smith 
8361321219cSEthan Coon 	lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
8371321219cSEthan Coon 	lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
83847c6ae99SBarry Smith 
83947c6ae99SBarry Smith 	cnt  = 0;
84047c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
84147c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
84247c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
843aa219208SBarry Smith 	      if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
84447c6ae99SBarry Smith 		cols[cnt++]  = k + nc*(slot + gnx*l + p);
84547c6ae99SBarry Smith 	      }
84647c6ae99SBarry Smith 	    }
84747c6ae99SBarry Smith 	  }
84847c6ae99SBarry Smith 	  rows[k]      = k + nc*(slot);
84947c6ae99SBarry Smith 	}
85047c6ae99SBarry Smith 	ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
85147c6ae99SBarry Smith       }
85247c6ae99SBarry Smith     }
85347c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
85447c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
85547c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
85647c6ae99SBarry Smith   }
85747c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
85847c6ae99SBarry Smith   PetscFunctionReturn(0);
85947c6ae99SBarry Smith }
86047c6ae99SBarry Smith 
86147c6ae99SBarry Smith #undef __FUNCT__
86294013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIAIJ_Fill"
86394013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J)
86447c6ae99SBarry Smith {
86547c6ae99SBarry Smith   PetscErrorCode         ierr;
86647c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
86747c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p;
86847c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
86947c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
87047c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
87147c6ae99SBarry Smith   MPI_Comm               comm;
87247c6ae99SBarry Smith   PetscScalar            *values;
8731321219cSEthan Coon   DMDABoundaryType       bx,by;
87447c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
875aa219208SBarry Smith   DMDAStencilType        st;
87647c6ae99SBarry Smith 
87747c6ae99SBarry Smith   PetscFunctionBegin;
87847c6ae99SBarry Smith   /*
87947c6ae99SBarry Smith          nc - number of components per grid point
88047c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
88147c6ae99SBarry Smith 
88247c6ae99SBarry Smith   */
8831321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
88447c6ae99SBarry Smith   col = 2*s + 1;
885aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
886aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
88747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
88847c6ae99SBarry Smith 
88947c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
8901411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
8911411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
89247c6ae99SBarry Smith 
89347c6ae99SBarry Smith   /* determine the matrix preallocation information */
89447c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
89547c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
89647c6ae99SBarry Smith 
8971321219cSEthan Coon     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
8981321219cSEthan Coon     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
89947c6ae99SBarry Smith 
90047c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
90147c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
90247c6ae99SBarry Smith 
9031321219cSEthan Coon       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9041321219cSEthan Coon       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
90547c6ae99SBarry Smith 
90647c6ae99SBarry Smith       for (k=0; k<nc; k++) {
90747c6ae99SBarry Smith         cnt  = 0;
90847c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
90947c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
91047c6ae99SBarry Smith             if (l || p) {
911aa219208SBarry Smith 	      if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
91247c6ae99SBarry Smith                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
91347c6ae99SBarry Smith 		  cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
91447c6ae99SBarry Smith 	      }
91547c6ae99SBarry Smith             } else {
91647c6ae99SBarry Smith 	      if (dfill) {
91747c6ae99SBarry Smith 		for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
91847c6ae99SBarry Smith 		  cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
91947c6ae99SBarry Smith 	      } else {
92047c6ae99SBarry Smith 		for (ifill_col=0; ifill_col<nc; ifill_col++)
92147c6ae99SBarry Smith 		  cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
92247c6ae99SBarry Smith 	      }
92347c6ae99SBarry Smith             }
92447c6ae99SBarry Smith 	  }
92547c6ae99SBarry Smith 	}
92647c6ae99SBarry Smith 	row = k + nc*(slot);
927784ac674SJed Brown         ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
92847c6ae99SBarry Smith       }
92947c6ae99SBarry Smith     }
93047c6ae99SBarry Smith   }
93147c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
93247c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
93347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
934784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
935784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
93647c6ae99SBarry Smith 
93747c6ae99SBarry Smith   /*
93847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
93947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
94047c6ae99SBarry Smith     PETSc ordering.
94147c6ae99SBarry Smith   */
94247c6ae99SBarry Smith   if (!dd->prealloc_only) {
94347c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
94447c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
94547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
94647c6ae99SBarry Smith 
9471321219cSEthan Coon       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
9481321219cSEthan Coon       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
94947c6ae99SBarry Smith 
95047c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
95147c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
95247c6ae99SBarry Smith 
9531321219cSEthan Coon 	lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
9541321219cSEthan Coon 	lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
95547c6ae99SBarry Smith 
95647c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
95747c6ae99SBarry Smith 	  cnt  = 0;
95847c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
95947c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
96047c6ae99SBarry Smith 	      if (l || p) {
961aa219208SBarry Smith 		if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
96247c6ae99SBarry Smith 		  for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
96347c6ae99SBarry Smith 		    cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
96447c6ae99SBarry Smith 		}
96547c6ae99SBarry Smith 	      } else {
96647c6ae99SBarry Smith 		if (dfill) {
96747c6ae99SBarry Smith 		  for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
96847c6ae99SBarry Smith 		    cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
96947c6ae99SBarry Smith 		} else {
97047c6ae99SBarry Smith 		  for (ifill_col=0; ifill_col<nc; ifill_col++)
97147c6ae99SBarry Smith 		    cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
97247c6ae99SBarry Smith 		}
97347c6ae99SBarry Smith 	      }
97447c6ae99SBarry Smith 	    }
97547c6ae99SBarry Smith 	  }
97647c6ae99SBarry Smith 	  row  = k + nc*(slot);
97747c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
97847c6ae99SBarry Smith 	}
97947c6ae99SBarry Smith       }
98047c6ae99SBarry Smith     }
98147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
98247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
98347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
98447c6ae99SBarry Smith   }
98547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
98647c6ae99SBarry Smith   PetscFunctionReturn(0);
98747c6ae99SBarry Smith }
98847c6ae99SBarry Smith 
98947c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
99047c6ae99SBarry Smith 
99147c6ae99SBarry Smith #undef __FUNCT__
99294013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIAIJ"
99394013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ(DM da,Mat J)
99447c6ae99SBarry Smith {
99547c6ae99SBarry Smith   PetscErrorCode         ierr;
99647c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
99747c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL;
99847c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
99947c6ae99SBarry Smith   MPI_Comm               comm;
100047c6ae99SBarry Smith   PetscScalar            *values;
10011321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
100247c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1003aa219208SBarry Smith   DMDAStencilType        st;
100447c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
100547c6ae99SBarry Smith 
100647c6ae99SBarry Smith   PetscFunctionBegin;
100747c6ae99SBarry Smith   /*
100847c6ae99SBarry Smith          nc - number of components per grid point
100947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
101047c6ae99SBarry Smith 
101147c6ae99SBarry Smith   */
10121321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
101347c6ae99SBarry Smith   col    = 2*s + 1;
101447c6ae99SBarry Smith 
1015aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1016aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
101747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
101847c6ae99SBarry Smith 
101947c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
10201411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
10211411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
102247c6ae99SBarry Smith 
102347c6ae99SBarry Smith   /* determine the matrix preallocation information */
102447c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
102547c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
10261321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10271321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
102847c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
10291321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10301321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
103147c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
10321321219cSEthan Coon 	kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
10331321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
103447c6ae99SBarry Smith 
103547c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
103647c6ae99SBarry Smith 
103747c6ae99SBarry Smith 	cnt  = 0;
103847c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
103947c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
104047c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
104147c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
1042aa219208SBarry Smith 		if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
104347c6ae99SBarry Smith 		  cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
104447c6ae99SBarry Smith 		}
104547c6ae99SBarry Smith 	      }
104647c6ae99SBarry Smith 	    }
104747c6ae99SBarry Smith 	  }
104847c6ae99SBarry Smith 	  rows[l] = l + nc*(slot);
104947c6ae99SBarry Smith 	}
1050784ac674SJed Brown 	ierr = MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
105147c6ae99SBarry Smith       }
105247c6ae99SBarry Smith     }
105347c6ae99SBarry Smith   }
105447c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
105547c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
105647c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
105747c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
1058784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1059784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
106047c6ae99SBarry Smith 
106147c6ae99SBarry Smith   /*
106247c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
106347c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
106447c6ae99SBarry Smith     PETSc ordering.
106547c6ae99SBarry Smith   */
106647c6ae99SBarry Smith   if (!dd->prealloc_only) {
106747c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
106847c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
106947c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
10701321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
10711321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
107247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
10731321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
10741321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
107547c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
10761321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
10771321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
107847c6ae99SBarry Smith 
107947c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
108047c6ae99SBarry Smith 
108147c6ae99SBarry Smith 	  cnt  = 0;
108247c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
108347c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
108447c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
108547c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
1086aa219208SBarry Smith 		  if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
108747c6ae99SBarry Smith 		    cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
108847c6ae99SBarry Smith 		  }
108947c6ae99SBarry Smith 		}
109047c6ae99SBarry Smith 	      }
109147c6ae99SBarry Smith 	    }
109247c6ae99SBarry Smith 	    rows[l]      = l + nc*(slot);
109347c6ae99SBarry Smith 	  }
109447c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
109547c6ae99SBarry Smith 	}
109647c6ae99SBarry Smith       }
109747c6ae99SBarry Smith     }
109847c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
109947c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
110047c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
110147c6ae99SBarry Smith   }
110247c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
110347c6ae99SBarry Smith   PetscFunctionReturn(0);
110447c6ae99SBarry Smith }
110547c6ae99SBarry Smith 
110647c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
110747c6ae99SBarry Smith 
110847c6ae99SBarry Smith #undef __FUNCT__
110994013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_1d_MPIAIJ"
111094013140SBarry Smith PetscErrorCode DMGetMatrix_DA_1d_MPIAIJ(DM da,Mat J)
111147c6ae99SBarry Smith {
111247c6ae99SBarry Smith   PetscErrorCode         ierr;
111347c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
111447c6ae99SBarry Smith   PetscInt               m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l;
111547c6ae99SBarry Smith   PetscInt               istart,iend;
111647c6ae99SBarry Smith   PetscScalar            *values;
11171321219cSEthan Coon   DMDABoundaryType       bx;
111847c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
111947c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
112047c6ae99SBarry Smith 
112147c6ae99SBarry Smith   PetscFunctionBegin;
112247c6ae99SBarry Smith   /*
112347c6ae99SBarry Smith          nc - number of components per grid point
112447c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
112547c6ae99SBarry Smith 
112647c6ae99SBarry Smith   */
11271321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);CHKERRQ(ierr);
112847c6ae99SBarry Smith   col    = 2*s + 1;
112947c6ae99SBarry Smith 
1130aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1131aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
113247c6ae99SBarry Smith 
113347c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
113447c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
113547c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
113647c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
113747c6ae99SBarry Smith 
11381411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
11391411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
1140784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1141784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
114247c6ae99SBarry Smith 
114347c6ae99SBarry Smith   /*
114447c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
114547c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
114647c6ae99SBarry Smith     PETSc ordering.
114747c6ae99SBarry Smith   */
114847c6ae99SBarry Smith   if (!dd->prealloc_only) {
114947c6ae99SBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
115047c6ae99SBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
115147c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
115247c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
115347c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
115447c6ae99SBarry Smith       slot   = i - gxs;
115547c6ae99SBarry Smith 
115647c6ae99SBarry Smith       cnt  = 0;
115747c6ae99SBarry Smith       for (l=0; l<nc; l++) {
115847c6ae99SBarry Smith 	for (i1=istart; i1<iend+1; i1++) {
115947c6ae99SBarry Smith 	  cols[cnt++] = l + nc*(slot + i1);
116047c6ae99SBarry Smith 	}
116147c6ae99SBarry Smith 	rows[l]      = l + nc*(slot);
116247c6ae99SBarry Smith       }
116347c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
116447c6ae99SBarry Smith     }
116547c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
116647c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
116747c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
116847c6ae99SBarry Smith   }
116947c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
117047c6ae99SBarry Smith   PetscFunctionReturn(0);
117147c6ae99SBarry Smith }
117247c6ae99SBarry Smith 
117347c6ae99SBarry Smith #undef __FUNCT__
117494013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPIBAIJ"
117594013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPIBAIJ(DM da,Mat J)
117647c6ae99SBarry Smith {
117747c6ae99SBarry Smith   PetscErrorCode         ierr;
117847c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
117947c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
118047c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
118147c6ae99SBarry Smith   MPI_Comm               comm;
118247c6ae99SBarry Smith   PetscScalar            *values;
11831321219cSEthan Coon   DMDABoundaryType       bx,by;
1184aa219208SBarry Smith   DMDAStencilType        st;
118547c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
118647c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
118747c6ae99SBarry Smith 
118847c6ae99SBarry Smith   PetscFunctionBegin;
118947c6ae99SBarry Smith   /*
119047c6ae99SBarry Smith      nc - number of components per grid point
119147c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
119247c6ae99SBarry Smith   */
11931321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
119447c6ae99SBarry Smith   col = 2*s + 1;
119547c6ae99SBarry Smith 
1196aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1197aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
119847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
119947c6ae99SBarry Smith 
120047c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
120147c6ae99SBarry Smith 
12021411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
12031411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
120447c6ae99SBarry Smith 
120547c6ae99SBarry Smith   /* determine the matrix preallocation information */
120647c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
120747c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
12081321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
12091321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
121047c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
12111321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
12121321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
121347c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
121447c6ae99SBarry Smith 
121547c6ae99SBarry Smith       /* Find block columns in block row */
121647c6ae99SBarry Smith       cnt  = 0;
121747c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
121847c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1219aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
122047c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
122147c6ae99SBarry Smith           }
122247c6ae99SBarry Smith         }
122347c6ae99SBarry Smith       }
1224784ac674SJed Brown       ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
122547c6ae99SBarry Smith     }
122647c6ae99SBarry Smith   }
122747c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
122847c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
122947c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
123047c6ae99SBarry Smith 
1231784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1232784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
123347c6ae99SBarry Smith 
123447c6ae99SBarry Smith   /*
123547c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
123647c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
123747c6ae99SBarry Smith     PETSc ordering.
123847c6ae99SBarry Smith   */
123947c6ae99SBarry Smith   if (!dd->prealloc_only) {
124047c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
124147c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
124247c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
12431321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
12441321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
124547c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
12461321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
12471321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
124847c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
124947c6ae99SBarry Smith 	cnt  = 0;
125047c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
125147c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1252aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
125347c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
125447c6ae99SBarry Smith             }
125547c6ae99SBarry Smith           }
125647c6ae99SBarry Smith         }
125747c6ae99SBarry Smith 	ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
125847c6ae99SBarry Smith       }
125947c6ae99SBarry Smith     }
126047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
126147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
126247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
126347c6ae99SBarry Smith   }
126447c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
126547c6ae99SBarry Smith   PetscFunctionReturn(0);
126647c6ae99SBarry Smith }
126747c6ae99SBarry Smith 
126847c6ae99SBarry Smith #undef __FUNCT__
126994013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIBAIJ"
127094013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIBAIJ(DM da,Mat J)
127147c6ae99SBarry Smith {
127247c6ae99SBarry Smith   PetscErrorCode         ierr;
127347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
127447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
127547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
127647c6ae99SBarry Smith   MPI_Comm               comm;
127747c6ae99SBarry Smith   PetscScalar            *values;
12781321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1279aa219208SBarry Smith   DMDAStencilType        st;
128047c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
128147c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
128247c6ae99SBarry Smith 
128347c6ae99SBarry Smith   PetscFunctionBegin;
128447c6ae99SBarry Smith   /*
128547c6ae99SBarry Smith          nc - number of components per grid point
128647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
128747c6ae99SBarry Smith 
128847c6ae99SBarry Smith   */
12891321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
129047c6ae99SBarry Smith   col    = 2*s + 1;
129147c6ae99SBarry Smith 
1292aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1293aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
129447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
129547c6ae99SBarry Smith 
129647c6ae99SBarry Smith   ierr  = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
129747c6ae99SBarry Smith 
12981411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
12991411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
130047c6ae99SBarry Smith 
130147c6ae99SBarry Smith   /* determine the matrix preallocation information */
130247c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
130347c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
13041321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
13051321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
130647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
13071321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
13081321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
130947c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
13101321219cSEthan Coon 	kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
13111321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
131247c6ae99SBarry Smith 
131347c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
131447c6ae99SBarry Smith 
131547c6ae99SBarry Smith 	/* Find block columns in block row */
131647c6ae99SBarry Smith 	cnt  = 0;
131747c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
131847c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
131947c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1320aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
132147c6ae99SBarry Smith 		cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
132247c6ae99SBarry Smith 	      }
132347c6ae99SBarry Smith 	    }
132447c6ae99SBarry Smith 	  }
132547c6ae99SBarry Smith 	}
1326784ac674SJed Brown 	ierr = MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);CHKERRQ(ierr);
132747c6ae99SBarry Smith       }
132847c6ae99SBarry Smith     }
132947c6ae99SBarry Smith   }
133047c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
133147c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
133247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
133347c6ae99SBarry Smith 
1334784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1335784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
133647c6ae99SBarry Smith 
133747c6ae99SBarry Smith   /*
133847c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
133947c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
134047c6ae99SBarry Smith     PETSc ordering.
134147c6ae99SBarry Smith   */
134247c6ae99SBarry Smith   if (!dd->prealloc_only) {
134347c6ae99SBarry Smith     ierr  = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
134447c6ae99SBarry Smith     ierr  = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
134547c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
13461321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
13471321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
134847c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
13491321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
13501321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
135147c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
13521321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
13531321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
135447c6ae99SBarry Smith 
135547c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
135647c6ae99SBarry Smith 
135747c6ae99SBarry Smith 	  cnt  = 0;
135847c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
135947c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
136047c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1361aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
136247c6ae99SBarry Smith                   cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
136347c6ae99SBarry Smith                 }
136447c6ae99SBarry Smith               }
136547c6ae99SBarry Smith             }
136647c6ae99SBarry Smith           }
136747c6ae99SBarry Smith 	  ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
136847c6ae99SBarry Smith 	}
136947c6ae99SBarry Smith       }
137047c6ae99SBarry Smith     }
137147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
137247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
137347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
137447c6ae99SBarry Smith   }
137547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
137647c6ae99SBarry Smith   PetscFunctionReturn(0);
137747c6ae99SBarry Smith }
137847c6ae99SBarry Smith 
137947c6ae99SBarry Smith #undef __FUNCT__
138047c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular"
138147c6ae99SBarry Smith /*
138247c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
138347c6ae99SBarry Smith   identify in the local ordering with periodic domain.
138447c6ae99SBarry Smith */
138547c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
138647c6ae99SBarry Smith {
138747c6ae99SBarry Smith   PetscErrorCode ierr;
138847c6ae99SBarry Smith   PetscInt       i,n;
138947c6ae99SBarry Smith 
139047c6ae99SBarry Smith   PetscFunctionBegin;
139147c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr);
139247c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr);
139347c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
139447c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
139547c6ae99SBarry Smith   }
139647c6ae99SBarry Smith   *cnt = n;
139747c6ae99SBarry Smith   PetscFunctionReturn(0);
139847c6ae99SBarry Smith }
139947c6ae99SBarry Smith 
140047c6ae99SBarry Smith #undef __FUNCT__
140194013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_2d_MPISBAIJ"
140294013140SBarry Smith PetscErrorCode DMGetMatrix_DA_2d_MPISBAIJ(DM da,Mat J)
140347c6ae99SBarry Smith {
140447c6ae99SBarry Smith   PetscErrorCode         ierr;
140547c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
140647c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
140747c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
140847c6ae99SBarry Smith   MPI_Comm               comm;
140947c6ae99SBarry Smith   PetscScalar            *values;
14101321219cSEthan Coon   DMDABoundaryType       bx,by;
1411aa219208SBarry Smith   DMDAStencilType        st;
141247c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
141347c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
141447c6ae99SBarry Smith 
141547c6ae99SBarry Smith   PetscFunctionBegin;
141647c6ae99SBarry Smith   /*
141747c6ae99SBarry Smith      nc - number of components per grid point
141847c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
141947c6ae99SBarry Smith   */
14201321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);CHKERRQ(ierr);
142147c6ae99SBarry Smith   col = 2*s + 1;
142247c6ae99SBarry Smith 
1423aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1424aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
142547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
142647c6ae99SBarry Smith 
142747c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
142847c6ae99SBarry Smith 
14291411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
14301411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
143147c6ae99SBarry Smith 
143247c6ae99SBarry Smith   /* determine the matrix preallocation information */
143347c6ae99SBarry Smith   ierr = MatPreallocateSymmetricInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
143447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
14351321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14361321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
143747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
14381321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14391321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
144047c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
144147c6ae99SBarry Smith 
144247c6ae99SBarry Smith       /* Find block columns in block row */
144347c6ae99SBarry Smith       cnt  = 0;
144447c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
144547c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1446aa219208SBarry Smith           if (st == DMDA_STENCIL_BOX || !ii || !jj) {
144747c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
144847c6ae99SBarry Smith           }
144947c6ae99SBarry Smith         }
145047c6ae99SBarry Smith       }
145147c6ae99SBarry Smith       ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
145247c6ae99SBarry Smith       ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
145347c6ae99SBarry Smith     }
145447c6ae99SBarry Smith   }
145547c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
145647c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
145747c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
145847c6ae99SBarry Smith 
1459784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1460784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
146147c6ae99SBarry Smith 
146247c6ae99SBarry Smith   /*
146347c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
146447c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
146547c6ae99SBarry Smith     PETSc ordering.
146647c6ae99SBarry Smith   */
146747c6ae99SBarry Smith   if (!dd->prealloc_only) {
146847c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
146947c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
147047c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
14711321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
14721321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
147347c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
14741321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
14751321219cSEthan Coon         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
147647c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
147747c6ae99SBarry Smith 
147847c6ae99SBarry Smith         /* Find block columns in block row */
147947c6ae99SBarry Smith         cnt  = 0;
148047c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
148147c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1482aa219208SBarry Smith             if (st == DMDA_STENCIL_BOX || !ii || !jj) {
148347c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
148447c6ae99SBarry Smith             }
148547c6ae99SBarry Smith           }
148647c6ae99SBarry Smith         }
148747c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
148847c6ae99SBarry Smith 	ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
148947c6ae99SBarry Smith       }
149047c6ae99SBarry Smith     }
149147c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
149247c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
149347c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
149447c6ae99SBarry Smith   }
149547c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
149647c6ae99SBarry Smith   PetscFunctionReturn(0);
149747c6ae99SBarry Smith }
149847c6ae99SBarry Smith 
149947c6ae99SBarry Smith #undef __FUNCT__
150094013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPISBAIJ"
150194013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPISBAIJ(DM da,Mat J)
150247c6ae99SBarry Smith {
150347c6ae99SBarry Smith   PetscErrorCode         ierr;
150447c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
150547c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
150647c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
150747c6ae99SBarry Smith   MPI_Comm               comm;
150847c6ae99SBarry Smith   PetscScalar            *values;
15091321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
1510aa219208SBarry Smith   DMDAStencilType        st;
151147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
151247c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
151347c6ae99SBarry Smith 
151447c6ae99SBarry Smith   PetscFunctionBegin;
151547c6ae99SBarry Smith   /*
151647c6ae99SBarry Smith      nc - number of components per grid point
151747c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
151847c6ae99SBarry Smith   */
15191321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
152047c6ae99SBarry Smith   col = 2*s + 1;
152147c6ae99SBarry Smith 
1522aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1523aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
152447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
152547c6ae99SBarry Smith 
152647c6ae99SBarry Smith   /* create the matrix */
152747c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
152847c6ae99SBarry Smith 
15291411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
15301411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
153147c6ae99SBarry Smith 
153247c6ae99SBarry Smith   /* determine the matrix preallocation information */
153347c6ae99SBarry Smith   ierr = MatPreallocateSymmetricInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
153447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
15351321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
15361321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
153747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
15381321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15391321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
154047c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
15411321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15421321219cSEthan Coon 	kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
154347c6ae99SBarry Smith 
154447c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
154547c6ae99SBarry Smith 
154647c6ae99SBarry Smith 	/* Find block columns in block row */
154747c6ae99SBarry Smith 	cnt  = 0;
154847c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
154947c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
155047c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1551aa219208SBarry Smith               if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
155247c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
155347c6ae99SBarry Smith               }
155447c6ae99SBarry Smith             }
155547c6ae99SBarry Smith           }
155647c6ae99SBarry Smith         }
155747c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
155847c6ae99SBarry Smith         ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
155947c6ae99SBarry Smith       }
156047c6ae99SBarry Smith     }
156147c6ae99SBarry Smith   }
156247c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
156347c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
156447c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
156547c6ae99SBarry Smith 
1566784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1567784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
156847c6ae99SBarry Smith 
156947c6ae99SBarry Smith   /*
157047c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
157147c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
157247c6ae99SBarry Smith     PETSc ordering.
157347c6ae99SBarry Smith   */
157447c6ae99SBarry Smith   if (!dd->prealloc_only) {
157547c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
157647c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
157747c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
15781321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
15791321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
158047c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
15811321219cSEthan Coon         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
15821321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
158347c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
15841321219cSEthan Coon           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
15851321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
158647c6ae99SBarry Smith 
158747c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
158847c6ae99SBarry Smith 
158947c6ae99SBarry Smith 	  cnt  = 0;
159047c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
159147c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
159247c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1593aa219208SBarry Smith                 if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
159447c6ae99SBarry Smith 		  cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
159547c6ae99SBarry Smith 		}
159647c6ae99SBarry Smith 	      }
159747c6ae99SBarry Smith 	    }
159847c6ae99SBarry Smith 	  }
159947c6ae99SBarry Smith           ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
160047c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
160147c6ae99SBarry Smith 	}
160247c6ae99SBarry Smith       }
160347c6ae99SBarry Smith     }
160447c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
160547c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160647c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160747c6ae99SBarry Smith   }
160847c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
160947c6ae99SBarry Smith   PetscFunctionReturn(0);
161047c6ae99SBarry Smith }
161147c6ae99SBarry Smith 
161247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
161347c6ae99SBarry Smith 
161447c6ae99SBarry Smith #undef __FUNCT__
161594013140SBarry Smith #define __FUNCT__ "DMGetMatrix_DA_3d_MPIAIJ_Fill"
161694013140SBarry Smith PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J)
161747c6ae99SBarry Smith {
161847c6ae99SBarry Smith   PetscErrorCode         ierr;
161947c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
162047c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz;
162147c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
162247c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
162347c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
162447c6ae99SBarry Smith   MPI_Comm               comm;
162547c6ae99SBarry Smith   PetscScalar            *values;
16261321219cSEthan Coon   DMDABoundaryType       bx,by,bz;
162747c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1628aa219208SBarry Smith   DMDAStencilType        st;
162947c6ae99SBarry Smith 
163047c6ae99SBarry Smith   PetscFunctionBegin;
163147c6ae99SBarry Smith   /*
163247c6ae99SBarry Smith          nc - number of components per grid point
163347c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
163447c6ae99SBarry Smith 
163547c6ae99SBarry Smith   */
16361321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);CHKERRQ(ierr);
163747c6ae99SBarry Smith   col    = 2*s + 1;
16381321219cSEthan Coon   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){
163947c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
164047c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
164147c6ae99SBarry Smith   }
16421321219cSEthan Coon   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){
164347c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
164447c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
164547c6ae99SBarry Smith   }
16461321219cSEthan Coon   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)){
164747c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
164847c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
164947c6ae99SBarry Smith   }
165047c6ae99SBarry Smith 
1651aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1652aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
165347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
165447c6ae99SBarry Smith 
165547c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
16561411c6eeSJed Brown   ierr = DMGetLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
16571411c6eeSJed Brown   ierr = DMGetLocalToGlobalMappingBlock(da,&ltogb);CHKERRQ(ierr);
165847c6ae99SBarry Smith 
165947c6ae99SBarry Smith   /* determine the matrix preallocation information */
166047c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
166147c6ae99SBarry Smith 
166247c6ae99SBarry Smith 
166347c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
16641321219cSEthan Coon     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
16651321219cSEthan Coon     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
166647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
16671321219cSEthan Coon       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
16681321219cSEthan Coon       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
166947c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
16701321219cSEthan Coon         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
16711321219cSEthan Coon         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
167247c6ae99SBarry Smith 
167347c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
167447c6ae99SBarry Smith 
167547c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
167647c6ae99SBarry Smith 	  cnt  = 0;
167747c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
167847c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
167947c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
168047c6ae99SBarry Smith 		if (ii || jj || kk) {
1681aa219208SBarry Smith 		  if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
168247c6ae99SBarry Smith 		    for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
168347c6ae99SBarry Smith 		      cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
168447c6ae99SBarry Smith 		  }
168547c6ae99SBarry Smith 		} else {
168647c6ae99SBarry Smith 		  if (dfill) {
168747c6ae99SBarry Smith 		    for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
168847c6ae99SBarry Smith 		      cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
168947c6ae99SBarry Smith 		  } else {
169047c6ae99SBarry Smith 		    for (ifill_col=0; ifill_col<nc; ifill_col++)
169147c6ae99SBarry Smith 		      cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
169247c6ae99SBarry Smith 		  }
169347c6ae99SBarry Smith 		}
169447c6ae99SBarry Smith 	      }
169547c6ae99SBarry Smith 	    }
169647c6ae99SBarry Smith 	  }
169747c6ae99SBarry Smith 	  row  = l + nc*(slot);
1698784ac674SJed Brown 	  ierr = MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);CHKERRQ(ierr);
169947c6ae99SBarry Smith 	}
170047c6ae99SBarry Smith       }
170147c6ae99SBarry Smith     }
170247c6ae99SBarry Smith   }
170347c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
170447c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
170547c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1706784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(J,ltog,ltog);CHKERRQ(ierr);
1707784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);CHKERRQ(ierr);
170847c6ae99SBarry Smith 
170947c6ae99SBarry Smith   /*
171047c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
171147c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
171247c6ae99SBarry Smith     PETSc ordering.
171347c6ae99SBarry Smith   */
171447c6ae99SBarry Smith   if (!dd->prealloc_only) {
171547c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
171647c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
171747c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
17181321219cSEthan Coon       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
17191321219cSEthan Coon       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
172047c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
17211321219cSEthan Coon 	jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
17221321219cSEthan Coon 	jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
172347c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
17241321219cSEthan Coon 	  kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
17251321219cSEthan Coon 	  kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
172647c6ae99SBarry Smith 
172747c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
172847c6ae99SBarry Smith 
172947c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
173047c6ae99SBarry Smith 	    cnt  = 0;
173147c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
173247c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
173347c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
173447c6ae99SBarry Smith 		  if (ii || jj || kk) {
1735aa219208SBarry Smith 		    if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
173647c6ae99SBarry Smith 		      for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
173747c6ae99SBarry Smith 			cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
173847c6ae99SBarry Smith 		    }
173947c6ae99SBarry Smith 		  } else {
174047c6ae99SBarry Smith 		    if (dfill) {
174147c6ae99SBarry Smith 		      for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
174247c6ae99SBarry Smith 			cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
174347c6ae99SBarry Smith 		    } else {
174447c6ae99SBarry Smith 		      for (ifill_col=0; ifill_col<nc; ifill_col++)
174547c6ae99SBarry Smith 			cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
174647c6ae99SBarry Smith 		    }
174747c6ae99SBarry Smith 		  }
174847c6ae99SBarry Smith 		}
174947c6ae99SBarry Smith 	      }
175047c6ae99SBarry Smith 	    }
175147c6ae99SBarry Smith 	    row  = l + nc*(slot);
175247c6ae99SBarry Smith 	    ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
175347c6ae99SBarry Smith 	  }
175447c6ae99SBarry Smith 	}
175547c6ae99SBarry Smith       }
175647c6ae99SBarry Smith     }
175747c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
175847c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
175947c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
176047c6ae99SBarry Smith   }
176147c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
176247c6ae99SBarry Smith   PetscFunctionReturn(0);
176347c6ae99SBarry Smith }
1764