xref: /petsc/src/dm/impls/da/fdda.c (revision 9a42bb27a39f0cdf3306a1e22d33cd9809484eaa)
147c6ae99SBarry Smith #define PETSCDM_DLL
247c6ae99SBarry Smith 
347c6ae99SBarry Smith #include "private/daimpl.h" /*I      "petscda.h"     I*/
447c6ae99SBarry Smith #include "petscmat.h"         /*I      "petscmat.h"    I*/
547c6ae99SBarry Smith #include "private/matimpl.h"
647c6ae99SBarry Smith 
7*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetColoring1d_MPIAIJ(DM,ISColoringType,ISColoring *);
8*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetColoring2d_MPIAIJ(DM,ISColoringType,ISColoring *);
9*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetColoring2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring *);
10*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetColoring3d_MPIAIJ(DM,ISColoringType,ISColoring *);
1147c6ae99SBarry Smith 
1247c6ae99SBarry Smith /*
1347c6ae99SBarry Smith    For ghost i that may be negative or greater than the upper bound this
1447c6ae99SBarry Smith   maps it into the 0:m-1 range using periodicity
1547c6ae99SBarry Smith */
1647c6ae99SBarry Smith #define SetInRange(i,m) ((i < 0) ? m+i:((i >= m) ? i-m:i))
1747c6ae99SBarry Smith 
1847c6ae99SBarry Smith #undef __FUNCT__
1947c6ae99SBarry Smith #define __FUNCT__ "DASetBlockFills_Private"
2047c6ae99SBarry Smith static PetscErrorCode DASetBlockFills_Private(PetscInt *dfill,PetscInt w,PetscInt **rfill)
2147c6ae99SBarry Smith {
2247c6ae99SBarry Smith   PetscErrorCode ierr;
2347c6ae99SBarry Smith   PetscInt       i,j,nz,*fill;
2447c6ae99SBarry Smith 
2547c6ae99SBarry Smith   PetscFunctionBegin;
2647c6ae99SBarry Smith   if (!dfill) PetscFunctionReturn(0);
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith   /* count number nonzeros */
2947c6ae99SBarry Smith   nz = 0;
3047c6ae99SBarry Smith   for (i=0; i<w; i++) {
3147c6ae99SBarry Smith     for (j=0; j<w; j++) {
3247c6ae99SBarry Smith       if (dfill[w*i+j]) nz++;
3347c6ae99SBarry Smith     }
3447c6ae99SBarry Smith   }
3547c6ae99SBarry Smith   ierr = PetscMalloc((nz + w + 1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
3647c6ae99SBarry Smith   /* construct modified CSR storage of nonzero structure */
3747c6ae99SBarry Smith   nz = w + 1;
3847c6ae99SBarry Smith   for (i=0; i<w; i++) {
3947c6ae99SBarry Smith     fill[i] = nz;
4047c6ae99SBarry Smith     for (j=0; j<w; j++) {
4147c6ae99SBarry Smith       if (dfill[w*i+j]) {
4247c6ae99SBarry Smith 	fill[nz] = j;
4347c6ae99SBarry Smith 	nz++;
4447c6ae99SBarry Smith       }
4547c6ae99SBarry Smith     }
4647c6ae99SBarry Smith   }
4747c6ae99SBarry Smith   fill[w] = nz;
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith   *rfill = fill;
5047c6ae99SBarry Smith   PetscFunctionReturn(0);
5147c6ae99SBarry Smith }
5247c6ae99SBarry Smith 
5347c6ae99SBarry Smith #undef __FUNCT__
5447c6ae99SBarry Smith #define __FUNCT__ "DASetMatPreallocateOnly"
5547c6ae99SBarry Smith /*@
5647c6ae99SBarry Smith     DASetMatPreallocateOnly - When DAGetMatrix() is called the matrix will be properly
5747c6ae99SBarry Smith        preallocated but the nonzero structure and zero values will not be set.
5847c6ae99SBarry Smith 
5947c6ae99SBarry Smith     Logically Collective on DA
6047c6ae99SBarry Smith 
6147c6ae99SBarry Smith     Input Parameter:
6247c6ae99SBarry Smith +   da - the distributed array
6347c6ae99SBarry Smith -   only - PETSC_TRUE if only want preallocation
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith 
6647c6ae99SBarry Smith     Level: developer
6747c6ae99SBarry Smith 
6847c6ae99SBarry Smith .seealso DAGetMatrix(), DASetGetMatrix(), DASetBlockSize(), DASetBlockFills()
6947c6ae99SBarry Smith 
7047c6ae99SBarry Smith @*/
71*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetMatPreallocateOnly(DM da,PetscBool  only)
7247c6ae99SBarry Smith {
7347c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
7447c6ae99SBarry Smith 
7547c6ae99SBarry Smith   PetscFunctionBegin;
7647c6ae99SBarry Smith   dd->prealloc_only = only;
7747c6ae99SBarry Smith   PetscFunctionReturn(0);
7847c6ae99SBarry Smith }
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith #undef __FUNCT__
8147c6ae99SBarry Smith #define __FUNCT__ "DASetBlockFills"
8247c6ae99SBarry Smith /*@
8347c6ae99SBarry Smith     DASetBlockFills - Sets the fill pattern in each block for a multi-component problem
8447c6ae99SBarry Smith     of the matrix returned by DAGetMatrix().
8547c6ae99SBarry Smith 
8647c6ae99SBarry Smith     Logically Collective on DA
8747c6ae99SBarry Smith 
8847c6ae99SBarry Smith     Input Parameter:
8947c6ae99SBarry Smith +   da - the distributed array
9047c6ae99SBarry Smith .   dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block)
9147c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
9247c6ae99SBarry Smith 
9347c6ae99SBarry Smith 
9447c6ae99SBarry Smith     Level: developer
9547c6ae99SBarry Smith 
9647c6ae99SBarry Smith     Notes: This only makes sense when you are doing multicomponent problems but using the
9747c6ae99SBarry Smith        MPIAIJ matrix format
9847c6ae99SBarry Smith 
9947c6ae99SBarry Smith            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
10047c6ae99SBarry Smith        representing coupling and 0 entries for missing coupling. For example
10147c6ae99SBarry Smith $             dfill[9] = {1, 0, 0,
10247c6ae99SBarry Smith $                         1, 1, 0,
10347c6ae99SBarry Smith $                         0, 1, 1}
10447c6ae99SBarry Smith        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
10547c6ae99SBarry Smith        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
10647c6ae99SBarry Smith        diagonal block).
10747c6ae99SBarry Smith 
10847c6ae99SBarry Smith      DASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
10947c6ae99SBarry Smith      can be represented in the dfill, ofill format
11047c6ae99SBarry Smith 
11147c6ae99SBarry Smith    Contributed by Glenn Hammond
11247c6ae99SBarry Smith 
11347c6ae99SBarry Smith .seealso DAGetMatrix(), DASetGetMatrix(), DASetMatPreallocateOnly()
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith @*/
116*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetBlockFills(DM da,PetscInt *dfill,PetscInt *ofill)
11747c6ae99SBarry Smith {
11847c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
11947c6ae99SBarry Smith   PetscErrorCode ierr;
12047c6ae99SBarry Smith 
12147c6ae99SBarry Smith   PetscFunctionBegin;
12247c6ae99SBarry Smith   ierr = DASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr);
12347c6ae99SBarry Smith   ierr = DASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr);
12447c6ae99SBarry Smith   PetscFunctionReturn(0);
12547c6ae99SBarry Smith }
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith #undef __FUNCT__
12947c6ae99SBarry Smith #define __FUNCT__ "DAGetColoring"
13047c6ae99SBarry Smith /*@
13147c6ae99SBarry Smith     DAGetColoring - Gets the coloring required for computing the Jacobian via
13247c6ae99SBarry Smith     finite differences on a function defined using a stencil on the DA.
13347c6ae99SBarry Smith 
13447c6ae99SBarry Smith     Collective on DA
13547c6ae99SBarry Smith 
13647c6ae99SBarry Smith     Input Parameter:
13747c6ae99SBarry Smith +   da - the distributed array
13847c6ae99SBarry Smith .   ctype - IS_COLORING_GLOBAL or IS_COLORING_GHOSTED
13947c6ae99SBarry Smith -   mtype - either MATAIJ or MATBAIJ
14047c6ae99SBarry Smith 
14147c6ae99SBarry Smith     Output Parameters:
14247c6ae99SBarry Smith .   coloring - matrix coloring for use in computing Jacobians (or PETSC_NULL if not needed)
14347c6ae99SBarry Smith 
14447c6ae99SBarry Smith     Level: advanced
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith     Notes: These compute the graph coloring of the graph of A^{T}A. The coloring used
14747c6ae99SBarry Smith    for efficient (parallel or thread based) triangular solves etc is NOT
14847c6ae99SBarry Smith    available.
14947c6ae99SBarry Smith 
15047c6ae99SBarry Smith         For BAIJ matrices this colors the graph for the blocks, not for the individual matrix elements;
15147c6ae99SBarry Smith     the same as MatGetColoring().
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), ISColoringType, ISColoring
15447c6ae99SBarry Smith 
15547c6ae99SBarry Smith @*/
156*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetColoring(DM da,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
15747c6ae99SBarry Smith {
15847c6ae99SBarry Smith   PetscErrorCode ierr;
15947c6ae99SBarry Smith   PetscInt       dim,m,n,p,nc;
16047c6ae99SBarry Smith   DAPeriodicType wrap;
16147c6ae99SBarry Smith   MPI_Comm       comm;
16247c6ae99SBarry Smith   PetscMPIInt    size;
16347c6ae99SBarry Smith   PetscBool      isBAIJ;
16447c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
16547c6ae99SBarry Smith 
16647c6ae99SBarry Smith   PetscFunctionBegin;
16747c6ae99SBarry Smith   /*
16847c6ae99SBarry Smith                                   m
16947c6ae99SBarry Smith           ------------------------------------------------------
17047c6ae99SBarry Smith          |                                                     |
17147c6ae99SBarry Smith          |                                                     |
17247c6ae99SBarry Smith          |               ----------------------                |
17347c6ae99SBarry Smith          |               |                    |                |
17447c6ae99SBarry Smith       n  |           yn  |                    |                |
17547c6ae99SBarry Smith          |               |                    |                |
17647c6ae99SBarry Smith          |               .---------------------                |
17747c6ae99SBarry Smith          |             (xs,ys)     xn                          |
17847c6ae99SBarry Smith          |            .                                        |
17947c6ae99SBarry Smith          |         (gxs,gys)                                   |
18047c6ae99SBarry Smith          |                                                     |
18147c6ae99SBarry Smith           -----------------------------------------------------
18247c6ae99SBarry Smith   */
18347c6ae99SBarry Smith 
18447c6ae99SBarry Smith   /*
18547c6ae99SBarry Smith          nc - number of components per grid point
18647c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
18747c6ae99SBarry Smith 
18847c6ae99SBarry Smith   */
18947c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&wrap,0);CHKERRQ(ierr);
19047c6ae99SBarry Smith 
19147c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
19247c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
19347c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED){
19447c6ae99SBarry Smith     if (size == 1) {
19547c6ae99SBarry Smith       ctype = IS_COLORING_GLOBAL;
19647c6ae99SBarry Smith     } else if (dim > 1){
19747c6ae99SBarry Smith       if ((m==1 && DAXPeriodic(wrap)) || (n==1 && DAYPeriodic(wrap)) || (p==1 && DAZPeriodic(wrap))){
19847c6ae99SBarry 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");
19947c6ae99SBarry Smith       }
20047c6ae99SBarry Smith     }
20147c6ae99SBarry Smith   }
20247c6ae99SBarry Smith 
20347c6ae99SBarry Smith   /* Tell the DA it has 1 degree of freedom per grid point so that the coloring for BAIJ
20447c6ae99SBarry Smith      matrices is for the blocks, not the individual matrix elements  */
20547c6ae99SBarry Smith   ierr = PetscStrcmp(mtype,MATBAIJ,&isBAIJ);CHKERRQ(ierr);
20647c6ae99SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);}
20747c6ae99SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);}
20847c6ae99SBarry Smith   if (isBAIJ) {
20947c6ae99SBarry Smith     dd->w = 1;
21047c6ae99SBarry Smith     dd->xs = dd->xs/nc;
21147c6ae99SBarry Smith     dd->xe = dd->xe/nc;
21247c6ae99SBarry Smith     dd->Xs = dd->Xs/nc;
21347c6ae99SBarry Smith     dd->Xe = dd->Xe/nc;
21447c6ae99SBarry Smith   }
21547c6ae99SBarry Smith 
21647c6ae99SBarry Smith   /*
21747c6ae99SBarry Smith      We do not provide a getcoloring function in the DA operations because
21847c6ae99SBarry Smith    the basic DA does not know about matrices. We think of DA as being more
21947c6ae99SBarry Smith    more low-level then matrices.
22047c6ae99SBarry Smith   */
22147c6ae99SBarry Smith   if (dim == 1) {
22247c6ae99SBarry Smith     ierr = DAGetColoring1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
22347c6ae99SBarry Smith   } else if (dim == 2) {
22447c6ae99SBarry Smith     ierr =  DAGetColoring2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
22547c6ae99SBarry Smith   } else if (dim == 3) {
22647c6ae99SBarry Smith     ierr =  DAGetColoring3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
22747c6ae99SBarry Smith   } else {
22847c6ae99SBarry Smith     SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim);
22947c6ae99SBarry Smith   }
23047c6ae99SBarry Smith   if (isBAIJ) {
23147c6ae99SBarry Smith     dd->w = nc;
23247c6ae99SBarry Smith     dd->xs = dd->xs*nc;
23347c6ae99SBarry Smith     dd->xe = dd->xe*nc;
23447c6ae99SBarry Smith     dd->Xs = dd->Xs*nc;
23547c6ae99SBarry Smith     dd->Xe = dd->Xe*nc;
23647c6ae99SBarry Smith   }
23747c6ae99SBarry Smith   PetscFunctionReturn(0);
23847c6ae99SBarry Smith }
23947c6ae99SBarry Smith 
24047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
24147c6ae99SBarry Smith 
24247c6ae99SBarry Smith #undef __FUNCT__
24347c6ae99SBarry Smith #define __FUNCT__ "DAGetColoring2d_MPIAIJ"
244*9a42bb27SBarry Smith PetscErrorCode DAGetColoring2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
24547c6ae99SBarry Smith {
24647c6ae99SBarry Smith   PetscErrorCode         ierr;
24747c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col;
24847c6ae99SBarry Smith   PetscInt               ncolors;
24947c6ae99SBarry Smith   MPI_Comm               comm;
25047c6ae99SBarry Smith   DAPeriodicType         wrap;
25147c6ae99SBarry Smith   DAStencilType          st;
25247c6ae99SBarry Smith   ISColoringValue        *colors;
25347c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
25447c6ae99SBarry Smith 
25547c6ae99SBarry Smith   PetscFunctionBegin;
25647c6ae99SBarry Smith   /*
25747c6ae99SBarry Smith          nc - number of components per grid point
25847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
25947c6ae99SBarry Smith 
26047c6ae99SBarry Smith   */
26147c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
26247c6ae99SBarry Smith   col    = 2*s + 1;
26347c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
26447c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
26547c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
26647c6ae99SBarry Smith 
26747c6ae99SBarry Smith   /* special case as taught to us by Paul Hovland */
26847c6ae99SBarry Smith   if (st == DA_STENCIL_STAR && s == 1) {
26947c6ae99SBarry Smith     ierr = DAGetColoring2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
27047c6ae99SBarry Smith   } else {
27147c6ae99SBarry Smith 
27247c6ae99SBarry Smith     if (DAXPeriodic(wrap) && (m % col)){
27347c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\
27447c6ae99SBarry Smith                  by 2*stencil_width + 1 (%d)\n", m, col);
27547c6ae99SBarry Smith     }
27647c6ae99SBarry Smith     if (DAYPeriodic(wrap) && (n % col)){
27747c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\
27847c6ae99SBarry Smith                  by 2*stencil_width + 1 (%d)\n", n, col);
27947c6ae99SBarry Smith     }
28047c6ae99SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
28147c6ae99SBarry Smith       if (!dd->localcoloring) {
28247c6ae99SBarry Smith 	ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
28347c6ae99SBarry Smith 	ii = 0;
28447c6ae99SBarry Smith 	for (j=ys; j<ys+ny; j++) {
28547c6ae99SBarry Smith 	  for (i=xs; i<xs+nx; i++) {
28647c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
28747c6ae99SBarry Smith 	      colors[ii++] = k + nc*((i % col) + col*(j % col));
28847c6ae99SBarry Smith 	    }
28947c6ae99SBarry Smith 	  }
29047c6ae99SBarry Smith 	}
29147c6ae99SBarry Smith         ncolors = nc + nc*(col-1 + col*(col-1));
29247c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
29347c6ae99SBarry Smith       }
29447c6ae99SBarry Smith       *coloring = dd->localcoloring;
29547c6ae99SBarry Smith     } else if (ctype == IS_COLORING_GHOSTED) {
29647c6ae99SBarry Smith       if (!dd->ghostedcoloring) {
29747c6ae99SBarry Smith 	ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
29847c6ae99SBarry Smith 	ii = 0;
29947c6ae99SBarry Smith 	for (j=gys; j<gys+gny; j++) {
30047c6ae99SBarry Smith 	  for (i=gxs; i<gxs+gnx; i++) {
30147c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
30247c6ae99SBarry Smith 	      /* the complicated stuff is to handle periodic boundaries */
30347c6ae99SBarry Smith 	      colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
30447c6ae99SBarry Smith 	    }
30547c6ae99SBarry Smith 	  }
30647c6ae99SBarry Smith 	}
30747c6ae99SBarry Smith         ncolors = nc + nc*(col - 1 + col*(col-1));
30847c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
30947c6ae99SBarry Smith         /* PetscIntView(ncolors,(PetscInt *)colors,0); */
31047c6ae99SBarry Smith 
31147c6ae99SBarry Smith 	ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
31247c6ae99SBarry Smith       }
31347c6ae99SBarry Smith       *coloring = dd->ghostedcoloring;
31447c6ae99SBarry Smith     } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
31547c6ae99SBarry Smith   }
31647c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
31747c6ae99SBarry Smith   PetscFunctionReturn(0);
31847c6ae99SBarry Smith }
31947c6ae99SBarry Smith 
32047c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
32147c6ae99SBarry Smith 
32247c6ae99SBarry Smith #undef __FUNCT__
32347c6ae99SBarry Smith #define __FUNCT__ "DAGetColoring3d_MPIAIJ"
324*9a42bb27SBarry Smith PetscErrorCode DAGetColoring3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
32547c6ae99SBarry Smith {
32647c6ae99SBarry Smith   PetscErrorCode  ierr;
32747c6ae99SBarry 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;
32847c6ae99SBarry Smith   PetscInt        ncolors;
32947c6ae99SBarry Smith   MPI_Comm        comm;
33047c6ae99SBarry Smith   DAPeriodicType  wrap;
33147c6ae99SBarry Smith   DAStencilType   st;
33247c6ae99SBarry Smith   ISColoringValue *colors;
33347c6ae99SBarry Smith   DM_DA           *dd = (DM_DA*)da->data;
33447c6ae99SBarry Smith 
33547c6ae99SBarry Smith   PetscFunctionBegin;
33647c6ae99SBarry Smith   /*
33747c6ae99SBarry Smith          nc - number of components per grid point
33847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
33947c6ae99SBarry Smith 
34047c6ae99SBarry Smith   */
34147c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&wrap,&st);CHKERRQ(ierr);
34247c6ae99SBarry Smith   col    = 2*s + 1;
34347c6ae99SBarry Smith   if (DAXPeriodic(wrap) && (m % col)){
34447c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
34547c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
34647c6ae99SBarry Smith   }
34747c6ae99SBarry Smith   if (DAYPeriodic(wrap) && (n % col)){
34847c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
34947c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
35047c6ae99SBarry Smith   }
35147c6ae99SBarry Smith   if (DAZPeriodic(wrap) && (p % col)){
35247c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
35347c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
35447c6ae99SBarry Smith   }
35547c6ae99SBarry Smith 
35647c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
35747c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
35847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
35947c6ae99SBarry Smith 
36047c6ae99SBarry Smith   /* create the coloring */
36147c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
36247c6ae99SBarry Smith     if (!dd->localcoloring) {
36347c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
36447c6ae99SBarry Smith       ii = 0;
36547c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
36647c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
36747c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
36847c6ae99SBarry Smith             for (l=0; l<nc; l++) {
36947c6ae99SBarry Smith               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
37047c6ae99SBarry Smith             }
37147c6ae99SBarry Smith           }
37247c6ae99SBarry Smith         }
37347c6ae99SBarry Smith       }
37447c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
37547c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr);
37647c6ae99SBarry Smith     }
37747c6ae99SBarry Smith     *coloring = dd->localcoloring;
37847c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
37947c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
38047c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
38147c6ae99SBarry Smith       ii = 0;
38247c6ae99SBarry Smith       for (k=gzs; k<gzs+gnz; k++) {
38347c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
38447c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
38547c6ae99SBarry Smith             for (l=0; l<nc; l++) {
38647c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
38747c6ae99SBarry Smith               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
38847c6ae99SBarry Smith             }
38947c6ae99SBarry Smith           }
39047c6ae99SBarry Smith         }
39147c6ae99SBarry Smith       }
39247c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
39347c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
39447c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
39547c6ae99SBarry Smith     }
39647c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
39747c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
39847c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
39947c6ae99SBarry Smith   PetscFunctionReturn(0);
40047c6ae99SBarry Smith }
40147c6ae99SBarry Smith 
40247c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
40347c6ae99SBarry Smith 
40447c6ae99SBarry Smith #undef __FUNCT__
40547c6ae99SBarry Smith #define __FUNCT__ "DAGetColoring1d_MPIAIJ"
406*9a42bb27SBarry Smith PetscErrorCode DAGetColoring1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
40747c6ae99SBarry Smith {
40847c6ae99SBarry Smith   PetscErrorCode  ierr;
40947c6ae99SBarry Smith   PetscInt        xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
41047c6ae99SBarry Smith   PetscInt        ncolors;
41147c6ae99SBarry Smith   MPI_Comm        comm;
41247c6ae99SBarry Smith   DAPeriodicType  wrap;
41347c6ae99SBarry Smith   ISColoringValue *colors;
41447c6ae99SBarry Smith   DM_DA           *dd = (DM_DA*)da->data;
41547c6ae99SBarry Smith 
41647c6ae99SBarry Smith   PetscFunctionBegin;
41747c6ae99SBarry Smith   /*
41847c6ae99SBarry Smith          nc - number of components per grid point
41947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
42047c6ae99SBarry Smith 
42147c6ae99SBarry Smith   */
42247c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&wrap,0);CHKERRQ(ierr);
42347c6ae99SBarry Smith   col    = 2*s + 1;
42447c6ae99SBarry Smith 
42547c6ae99SBarry Smith   if (DAXPeriodic(wrap) && (m % col)) {
42647c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points is divisible\n\
42747c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
42847c6ae99SBarry Smith   }
42947c6ae99SBarry Smith 
43047c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
43147c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
43247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
43347c6ae99SBarry Smith 
43447c6ae99SBarry Smith   /* create the coloring */
43547c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
43647c6ae99SBarry Smith     if (!dd->localcoloring) {
43747c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
43847c6ae99SBarry Smith       i1 = 0;
43947c6ae99SBarry Smith       for (i=xs; i<xs+nx; i++) {
44047c6ae99SBarry Smith         for (l=0; l<nc; l++) {
44147c6ae99SBarry Smith           colors[i1++] = l + nc*(i % col);
44247c6ae99SBarry Smith         }
44347c6ae99SBarry Smith       }
44447c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
44547c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr);
44647c6ae99SBarry Smith     }
44747c6ae99SBarry Smith     *coloring = dd->localcoloring;
44847c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
44947c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
45047c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
45147c6ae99SBarry Smith       i1 = 0;
45247c6ae99SBarry Smith       for (i=gxs; i<gxs+gnx; i++) {
45347c6ae99SBarry Smith         for (l=0; l<nc; l++) {
45447c6ae99SBarry Smith           /* the complicated stuff is to handle periodic boundaries */
45547c6ae99SBarry Smith           colors[i1++] = l + nc*(SetInRange(i,m) % col);
45647c6ae99SBarry Smith         }
45747c6ae99SBarry Smith       }
45847c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
45947c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
46047c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
46147c6ae99SBarry Smith     }
46247c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
46347c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
46447c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
46547c6ae99SBarry Smith   PetscFunctionReturn(0);
46647c6ae99SBarry Smith }
46747c6ae99SBarry Smith 
46847c6ae99SBarry Smith #undef __FUNCT__
46947c6ae99SBarry Smith #define __FUNCT__ "DAGetColoring2d_5pt_MPIAIJ"
470*9a42bb27SBarry Smith PetscErrorCode DAGetColoring2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
47147c6ae99SBarry Smith {
47247c6ae99SBarry Smith   PetscErrorCode  ierr;
47347c6ae99SBarry Smith   PetscInt        xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
47447c6ae99SBarry Smith   PetscInt        ncolors;
47547c6ae99SBarry Smith   MPI_Comm        comm;
47647c6ae99SBarry Smith   DAPeriodicType  wrap;
47747c6ae99SBarry Smith   ISColoringValue *colors;
47847c6ae99SBarry Smith   DM_DA           *dd = (DM_DA*)da->data;
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith   PetscFunctionBegin;
48147c6ae99SBarry Smith   /*
48247c6ae99SBarry Smith          nc - number of components per grid point
48347c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
48447c6ae99SBarry Smith 
48547c6ae99SBarry Smith   */
48647c6ae99SBarry Smith   ierr   = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,0);CHKERRQ(ierr);
48747c6ae99SBarry Smith   ierr   = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
48847c6ae99SBarry Smith   ierr   = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
48947c6ae99SBarry Smith   ierr   = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
49047c6ae99SBarry Smith 
49147c6ae99SBarry Smith   if (DAXPeriodic(wrap) && (m % 5)){
49247c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
49347c6ae99SBarry Smith                  by 5\n");
49447c6ae99SBarry Smith   }
49547c6ae99SBarry Smith   if (DAYPeriodic(wrap) && (n % 5)){
49647c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
49747c6ae99SBarry Smith                  by 5\n");
49847c6ae99SBarry Smith   }
49947c6ae99SBarry Smith 
50047c6ae99SBarry Smith   /* create the coloring */
50147c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
50247c6ae99SBarry Smith     if (!dd->localcoloring) {
50347c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
50447c6ae99SBarry Smith       ii = 0;
50547c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
50647c6ae99SBarry Smith 	for (i=xs; i<xs+nx; i++) {
50747c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
50847c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*j+i) % 5);
50947c6ae99SBarry Smith 	  }
51047c6ae99SBarry Smith 	}
51147c6ae99SBarry Smith       }
51247c6ae99SBarry Smith       ncolors = 5*nc;
51347c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
51447c6ae99SBarry Smith     }
51547c6ae99SBarry Smith     *coloring = dd->localcoloring;
51647c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
51747c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
51847c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
51947c6ae99SBarry Smith       ii = 0;
52047c6ae99SBarry Smith       for (j=gys; j<gys+gny; j++) {
52147c6ae99SBarry Smith 	for (i=gxs; i<gxs+gnx; i++) {
52247c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
52347c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
52447c6ae99SBarry Smith 	  }
52547c6ae99SBarry Smith 	}
52647c6ae99SBarry Smith       }
52747c6ae99SBarry Smith       ncolors = 5*nc;
52847c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
52947c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
53047c6ae99SBarry Smith     }
53147c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
53247c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
53347c6ae99SBarry Smith   PetscFunctionReturn(0);
53447c6ae99SBarry Smith }
53547c6ae99SBarry Smith 
53647c6ae99SBarry Smith /* =========================================================================== */
537*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetMatrix1d_MPIAIJ(DM,Mat);
538*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetMatrix2d_MPIAIJ(DM,Mat);
539*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetMatrix2d_MPIAIJ_Fill(DM,Mat);
540*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetMatrix3d_MPIAIJ(DM,Mat);
541*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetMatrix3d_MPIAIJ_Fill(DM,Mat);
542*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetMatrix2d_MPIBAIJ(DM,Mat);
543*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetMatrix3d_MPIBAIJ(DM,Mat);
544*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetMatrix2d_MPISBAIJ(DM,Mat);
545*9a42bb27SBarry Smith EXTERN PetscErrorCode DAGetMatrix3d_MPISBAIJ(DM,Mat);
54647c6ae99SBarry Smith 
54747c6ae99SBarry Smith #undef __FUNCT__
54847c6ae99SBarry Smith #define __FUNCT__ "MatSetDA"
54947c6ae99SBarry Smith /*@
55047c6ae99SBarry Smith    MatSetDA - Sets the DA that is to be used by the HYPRE_StructMatrix PETSc matrix
55147c6ae99SBarry Smith 
55247c6ae99SBarry Smith    Logically Collective on Mat
55347c6ae99SBarry Smith 
55447c6ae99SBarry Smith    Input Parameters:
55547c6ae99SBarry Smith +  mat - the matrix
55647c6ae99SBarry Smith -  da - the da
55747c6ae99SBarry Smith 
55847c6ae99SBarry Smith    Level: intermediate
55947c6ae99SBarry Smith 
56047c6ae99SBarry Smith @*/
561*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT MatSetDA(Mat mat,DM da)
56247c6ae99SBarry Smith {
56347c6ae99SBarry Smith   PetscErrorCode ierr;
56447c6ae99SBarry Smith 
56547c6ae99SBarry Smith   PetscFunctionBegin;
56647c6ae99SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
56747c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
568*9a42bb27SBarry Smith   ierr = PetscTryMethod(mat,"MatSetDA_C",(Mat,DM),(mat,da));CHKERRQ(ierr);
56947c6ae99SBarry Smith   PetscFunctionReturn(0);
57047c6ae99SBarry Smith }
57147c6ae99SBarry Smith 
57247c6ae99SBarry Smith EXTERN_C_BEGIN
57347c6ae99SBarry Smith #undef __FUNCT__
57447c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA"
57547c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT MatView_MPI_DA(Mat A,PetscViewer viewer)
57647c6ae99SBarry Smith {
577*9a42bb27SBarry Smith   DM             da;
57847c6ae99SBarry Smith   PetscErrorCode ierr;
57947c6ae99SBarry Smith   const char     *prefix;
58047c6ae99SBarry Smith   Mat            Anatural;
58147c6ae99SBarry Smith   AO             ao;
58247c6ae99SBarry Smith   PetscInt       rstart,rend,*petsc,i;
58347c6ae99SBarry Smith   IS             is;
58447c6ae99SBarry Smith   MPI_Comm       comm;
58547c6ae99SBarry Smith 
58647c6ae99SBarry Smith   PetscFunctionBegin;
58747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
58847c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"DA",(PetscObject*)&da);CHKERRQ(ierr);
58947c6ae99SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DA");
59047c6ae99SBarry Smith 
59147c6ae99SBarry Smith   ierr = DAGetAO(da,&ao);CHKERRQ(ierr);
59247c6ae99SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
59347c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr);
59447c6ae99SBarry Smith   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
59547c6ae99SBarry Smith   ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr);
59647c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
59747c6ae99SBarry Smith 
59847c6ae99SBarry Smith   /* call viewer on natural ordering */
59947c6ae99SBarry Smith   ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr);
60047c6ae99SBarry Smith   ierr = ISDestroy(is);CHKERRQ(ierr);
60147c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr);
60247c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr);
60347c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr);
60447c6ae99SBarry Smith   ierr = MatView(Anatural,viewer);CHKERRQ(ierr);
60547c6ae99SBarry Smith   ierr = MatDestroy(Anatural);CHKERRQ(ierr);
60647c6ae99SBarry Smith   PetscFunctionReturn(0);
60747c6ae99SBarry Smith }
60847c6ae99SBarry Smith EXTERN_C_END
60947c6ae99SBarry Smith 
61047c6ae99SBarry Smith EXTERN_C_BEGIN
61147c6ae99SBarry Smith #undef __FUNCT__
61247c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA"
61347c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT MatLoad_MPI_DA(Mat A,PetscViewer viewer)
61447c6ae99SBarry Smith {
615*9a42bb27SBarry Smith   DM             da;
61647c6ae99SBarry Smith   PetscErrorCode ierr;
61747c6ae99SBarry Smith   Mat            Anatural,Aapp;
61847c6ae99SBarry Smith   AO             ao;
61947c6ae99SBarry Smith   PetscInt       rstart,rend,*app,i;
62047c6ae99SBarry Smith   IS             is;
62147c6ae99SBarry Smith   MPI_Comm       comm;
62247c6ae99SBarry Smith 
62347c6ae99SBarry Smith   PetscFunctionBegin;
62447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
62547c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"DA",(PetscObject*)&da);CHKERRQ(ierr);
62647c6ae99SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DA");
62747c6ae99SBarry Smith 
62847c6ae99SBarry Smith   /* Load the matrix in natural ordering */
62947c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,&Anatural);CHKERRQ(ierr);
63047c6ae99SBarry Smith   ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr);
63147c6ae99SBarry Smith   ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
63247c6ae99SBarry Smith   ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr);
63347c6ae99SBarry Smith 
63447c6ae99SBarry Smith   /* Map natural ordering to application ordering and create IS */
63547c6ae99SBarry Smith   ierr = DAGetAO(da,&ao);CHKERRQ(ierr);
63647c6ae99SBarry Smith   ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr);
63747c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr);
63847c6ae99SBarry Smith   for (i=rstart; i<rend; i++) app[i-rstart] = i;
63947c6ae99SBarry Smith   ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr);
64047c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
64147c6ae99SBarry Smith 
64247c6ae99SBarry Smith   /* Do permutation and replace header */
64347c6ae99SBarry Smith   ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr);
64447c6ae99SBarry Smith   ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr);
64547c6ae99SBarry Smith   ierr = ISDestroy(is);CHKERRQ(ierr);
64647c6ae99SBarry Smith   ierr = MatDestroy(Anatural);CHKERRQ(ierr);
64747c6ae99SBarry Smith   PetscFunctionReturn(0);
64847c6ae99SBarry Smith }
64947c6ae99SBarry Smith EXTERN_C_END
65047c6ae99SBarry Smith 
65147c6ae99SBarry Smith 
65247c6ae99SBarry Smith #undef __FUNCT__
65347c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix"
65447c6ae99SBarry Smith /*@C
65547c6ae99SBarry Smith     DAGetMatrix - Creates a matrix with the correct parallel layout and nonzero structure required for
65647c6ae99SBarry Smith       computing the Jacobian on a function defined using the stencil set in the DA.
65747c6ae99SBarry Smith 
65847c6ae99SBarry Smith     Collective on DA
65947c6ae99SBarry Smith 
66047c6ae99SBarry Smith     Input Parameter:
66147c6ae99SBarry Smith +   da - the distributed array
66247c6ae99SBarry Smith -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ,
66347c6ae99SBarry Smith             or any type which inherits from one of these (such as MATAIJ, MATBAIJ, MATSBAIJ)
66447c6ae99SBarry Smith 
66547c6ae99SBarry Smith     Output Parameters:
66647c6ae99SBarry Smith .   J  - matrix with the correct nonzero structure
66747c6ae99SBarry Smith         (obviously without the correct Jacobian values)
66847c6ae99SBarry Smith 
66947c6ae99SBarry Smith     Level: advanced
67047c6ae99SBarry Smith 
67147c6ae99SBarry Smith     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
67247c6ae99SBarry Smith        do not need to do it yourself.
67347c6ae99SBarry Smith 
67447c6ae99SBarry Smith        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
67547c6ae99SBarry Smith        the nonzero pattern call DASetMatPreallocateOnly()
67647c6ae99SBarry Smith 
67747c6ae99SBarry Smith        When you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
67847c6ae99SBarry Smith        internally by PETSc.
67947c6ae99SBarry Smith 
68047c6ae99SBarry Smith        In general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
68147c6ae99SBarry Smith        the indices for the global numbering for DAs which is complicated.
68247c6ae99SBarry Smith 
68347c6ae99SBarry Smith .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), DASetBlockFills(), DASetMatPreallocateOnly()
68447c6ae99SBarry Smith 
68547c6ae99SBarry Smith @*/
686*9a42bb27SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetMatrix(DM da, const MatType mtype,Mat *J)
68747c6ae99SBarry Smith {
68847c6ae99SBarry Smith   PetscErrorCode ierr;
68947c6ae99SBarry Smith   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
69047c6ae99SBarry Smith   Mat            A;
69147c6ae99SBarry Smith   MPI_Comm       comm;
69247c6ae99SBarry Smith   const MatType  Atype;
69347c6ae99SBarry Smith   void           (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL;
69447c6ae99SBarry Smith   MatType        ttype[256];
69547c6ae99SBarry Smith   PetscBool      flg;
69647c6ae99SBarry Smith   PetscMPIInt    size;
69747c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
69847c6ae99SBarry Smith 
69947c6ae99SBarry Smith   PetscFunctionBegin;
70047c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
70147c6ae99SBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
70247c6ae99SBarry Smith #endif
70347c6ae99SBarry Smith   ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr);
70447c6ae99SBarry Smith   ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DA options","Mat");CHKERRQ(ierr);
70547c6ae99SBarry Smith   ierr = PetscOptionsList("-da_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr);
70647c6ae99SBarry Smith   ierr = PetscOptionsEnd();
70747c6ae99SBarry Smith 
70847c6ae99SBarry Smith   /*
70947c6ae99SBarry Smith                                   m
71047c6ae99SBarry Smith           ------------------------------------------------------
71147c6ae99SBarry Smith          |                                                     |
71247c6ae99SBarry Smith          |                                                     |
71347c6ae99SBarry Smith          |               ----------------------                |
71447c6ae99SBarry Smith          |               |                    |                |
71547c6ae99SBarry Smith       n  |           ny  |                    |                |
71647c6ae99SBarry Smith          |               |                    |                |
71747c6ae99SBarry Smith          |               .---------------------                |
71847c6ae99SBarry Smith          |             (xs,ys)     nx                          |
71947c6ae99SBarry Smith          |            .                                        |
72047c6ae99SBarry Smith          |         (gxs,gys)                                   |
72147c6ae99SBarry Smith          |                                                     |
72247c6ae99SBarry Smith           -----------------------------------------------------
72347c6ae99SBarry Smith   */
72447c6ae99SBarry Smith 
72547c6ae99SBarry Smith   /*
72647c6ae99SBarry Smith          nc - number of components per grid point
72747c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
72847c6ae99SBarry Smith 
72947c6ae99SBarry Smith   */
73047c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0);CHKERRQ(ierr);
73147c6ae99SBarry Smith   ierr = DAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
73247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
73347c6ae99SBarry Smith   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
73447c6ae99SBarry Smith   ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr);
73547c6ae99SBarry Smith   ierr = MatSetType(A,(const MatType)ttype);CHKERRQ(ierr);
73647c6ae99SBarry Smith   ierr = MatSetDA(A,da);CHKERRQ(ierr);
73747c6ae99SBarry Smith   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
73847c6ae99SBarry Smith   ierr = MatGetType(A,&Atype);CHKERRQ(ierr);
73947c6ae99SBarry Smith   /*
74047c6ae99SBarry Smith      We do not provide a getmatrix function in the DA operations because
74147c6ae99SBarry Smith    the basic DA does not know about matrices. We think of DA as being more
74247c6ae99SBarry Smith    more low-level than matrices. This is kind of cheating but, cause sometimes
74347c6ae99SBarry Smith    we think of DA has higher level than matrices.
74447c6ae99SBarry Smith 
74547c6ae99SBarry Smith      We could switch based on Atype (or mtype), but we do not since the
74647c6ae99SBarry Smith    specialized setting routines depend only the particular preallocation
74747c6ae99SBarry Smith    details of the matrix, not the type itself.
74847c6ae99SBarry Smith   */
74947c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
75047c6ae99SBarry Smith   if (!aij) {
75147c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
75247c6ae99SBarry Smith   }
75347c6ae99SBarry Smith   if (!aij) {
75447c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
75547c6ae99SBarry Smith     if (!baij) {
75647c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
75747c6ae99SBarry Smith     }
75847c6ae99SBarry Smith     if (!baij){
75947c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
76047c6ae99SBarry Smith       if (!sbaij) {
76147c6ae99SBarry Smith         ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
76247c6ae99SBarry Smith       }
76347c6ae99SBarry Smith     }
76447c6ae99SBarry Smith   }
76547c6ae99SBarry Smith   if (aij) {
76647c6ae99SBarry Smith     if (dim == 1) {
76747c6ae99SBarry Smith       ierr = DAGetMatrix1d_MPIAIJ(da,A);CHKERRQ(ierr);
76847c6ae99SBarry Smith     } else if (dim == 2) {
76947c6ae99SBarry Smith       if (dd->ofill) {
77047c6ae99SBarry Smith         ierr = DAGetMatrix2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
77147c6ae99SBarry Smith       } else {
77247c6ae99SBarry Smith         ierr = DAGetMatrix2d_MPIAIJ(da,A);CHKERRQ(ierr);
77347c6ae99SBarry Smith       }
77447c6ae99SBarry Smith     } else if (dim == 3) {
77547c6ae99SBarry Smith       if (dd->ofill) {
77647c6ae99SBarry Smith         ierr = DAGetMatrix3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
77747c6ae99SBarry Smith       } else {
77847c6ae99SBarry Smith         ierr = DAGetMatrix3d_MPIAIJ(da,A);CHKERRQ(ierr);
77947c6ae99SBarry Smith       }
78047c6ae99SBarry Smith     }
78147c6ae99SBarry Smith   } else if (baij) {
78247c6ae99SBarry Smith     if (dim == 2) {
78347c6ae99SBarry Smith       ierr = DAGetMatrix2d_MPIBAIJ(da,A);CHKERRQ(ierr);
78447c6ae99SBarry Smith     } else if (dim == 3) {
78547c6ae99SBarry Smith       ierr = DAGetMatrix3d_MPIBAIJ(da,A);CHKERRQ(ierr);
78647c6ae99SBarry Smith     } else {
78747c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \
78847c6ae99SBarry Smith 	       "Send mail to petsc-maint@mcs.anl.gov for code",Atype,dim);
78947c6ae99SBarry Smith     }
79047c6ae99SBarry Smith   } else if (sbaij) {
79147c6ae99SBarry Smith     if (dim == 2) {
79247c6ae99SBarry Smith       ierr = DAGetMatrix2d_MPISBAIJ(da,A);CHKERRQ(ierr);
79347c6ae99SBarry Smith     } else if (dim == 3) {
79447c6ae99SBarry Smith       ierr = DAGetMatrix3d_MPISBAIJ(da,A);CHKERRQ(ierr);
79547c6ae99SBarry Smith     } else {
79647c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \
79747c6ae99SBarry Smith 	       "Send mail to petsc-maint@mcs.anl.gov for code",Atype,dim);
79847c6ae99SBarry Smith     }
79947c6ae99SBarry Smith   }
80047c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
80147c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
80247c6ae99SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"DA",(PetscObject)da);CHKERRQ(ierr);
80347c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
80447c6ae99SBarry Smith   if (size > 1) {
80547c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
80647c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);CHKERRQ(ierr);
80747c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);CHKERRQ(ierr);
80847c6ae99SBarry Smith   }
80947c6ae99SBarry Smith   *J = A;
81047c6ae99SBarry Smith   PetscFunctionReturn(0);
81147c6ae99SBarry Smith }
81247c6ae99SBarry Smith 
81347c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
81447c6ae99SBarry Smith #undef __FUNCT__
81547c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix2d_MPIAIJ"
816*9a42bb27SBarry Smith PetscErrorCode DAGetMatrix2d_MPIAIJ(DM da,Mat J)
81747c6ae99SBarry Smith {
81847c6ae99SBarry Smith   PetscErrorCode         ierr;
81947c6ae99SBarry 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;
82047c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
82147c6ae99SBarry Smith   MPI_Comm               comm;
82247c6ae99SBarry Smith   PetscScalar            *values;
82347c6ae99SBarry Smith   DAPeriodicType         wrap;
82447c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
82547c6ae99SBarry Smith   DAStencilType          st;
82647c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
82747c6ae99SBarry Smith 
82847c6ae99SBarry Smith   PetscFunctionBegin;
82947c6ae99SBarry Smith   /*
83047c6ae99SBarry Smith          nc - number of components per grid point
83147c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
83247c6ae99SBarry Smith 
83347c6ae99SBarry Smith   */
83447c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
83547c6ae99SBarry Smith   col = 2*s + 1;
83647c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
83747c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
83847c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
83947c6ae99SBarry Smith 
84047c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
84147c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
84247c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
84347c6ae99SBarry Smith 
84447c6ae99SBarry Smith   /* determine the matrix preallocation information */
84547c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
84647c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
84747c6ae99SBarry Smith 
84847c6ae99SBarry Smith     pstart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
84947c6ae99SBarry Smith     pend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
85047c6ae99SBarry Smith 
85147c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
85247c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
85347c6ae99SBarry Smith 
85447c6ae99SBarry Smith       lstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
85547c6ae99SBarry Smith       lend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
85647c6ae99SBarry Smith 
85747c6ae99SBarry Smith       cnt  = 0;
85847c6ae99SBarry Smith       for (k=0; k<nc; k++) {
85947c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
86047c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
86147c6ae99SBarry Smith 	    if ((st == DA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
86247c6ae99SBarry Smith 	      cols[cnt++]  = k + nc*(slot + gnx*l + p);
86347c6ae99SBarry Smith 	    }
86447c6ae99SBarry Smith 	  }
86547c6ae99SBarry Smith 	}
86647c6ae99SBarry Smith 	rows[k] = k + nc*(slot);
86747c6ae99SBarry Smith       }
86847c6ae99SBarry Smith       ierr = MatPreallocateSetLocal(ltog,nc,rows,cnt,cols,dnz,onz);CHKERRQ(ierr);
86947c6ae99SBarry Smith     }
87047c6ae99SBarry Smith   }
87147c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
87247c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
87347c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
87447c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
87547c6ae99SBarry Smith 
87647c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
87747c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
87847c6ae99SBarry Smith 
87947c6ae99SBarry Smith   /*
88047c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
88147c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
88247c6ae99SBarry Smith     PETSc ordering.
88347c6ae99SBarry Smith   */
88447c6ae99SBarry Smith   if (!dd->prealloc_only) {
88547c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
88647c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
88747c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
88847c6ae99SBarry Smith 
88947c6ae99SBarry Smith       pstart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
89047c6ae99SBarry Smith       pend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
89147c6ae99SBarry Smith 
89247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
89347c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
89447c6ae99SBarry Smith 
89547c6ae99SBarry Smith 	lstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
89647c6ae99SBarry Smith 	lend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
89747c6ae99SBarry Smith 
89847c6ae99SBarry Smith 	cnt  = 0;
89947c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
90047c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
90147c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
90247c6ae99SBarry Smith 	      if ((st == DA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
90347c6ae99SBarry Smith 		cols[cnt++]  = k + nc*(slot + gnx*l + p);
90447c6ae99SBarry Smith 	      }
90547c6ae99SBarry Smith 	    }
90647c6ae99SBarry Smith 	  }
90747c6ae99SBarry Smith 	  rows[k]      = k + nc*(slot);
90847c6ae99SBarry Smith 	}
90947c6ae99SBarry Smith 	ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
91047c6ae99SBarry Smith       }
91147c6ae99SBarry Smith     }
91247c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
91347c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91447c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91547c6ae99SBarry Smith   }
91647c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
91747c6ae99SBarry Smith   PetscFunctionReturn(0);
91847c6ae99SBarry Smith }
91947c6ae99SBarry Smith 
92047c6ae99SBarry Smith #undef __FUNCT__
92147c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix2d_MPIAIJ_Fill"
922*9a42bb27SBarry Smith PetscErrorCode DAGetMatrix2d_MPIAIJ_Fill(DM da,Mat J)
92347c6ae99SBarry Smith {
92447c6ae99SBarry Smith   PetscErrorCode         ierr;
92547c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
92647c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p;
92747c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
92847c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
92947c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
93047c6ae99SBarry Smith   MPI_Comm               comm;
93147c6ae99SBarry Smith   PetscScalar            *values;
93247c6ae99SBarry Smith   DAPeriodicType         wrap;
93347c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
93447c6ae99SBarry Smith   DAStencilType          st;
93547c6ae99SBarry Smith 
93647c6ae99SBarry Smith   PetscFunctionBegin;
93747c6ae99SBarry Smith   /*
93847c6ae99SBarry Smith          nc - number of components per grid point
93947c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
94047c6ae99SBarry Smith 
94147c6ae99SBarry Smith   */
94247c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
94347c6ae99SBarry Smith   col = 2*s + 1;
94447c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
94547c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
94647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
94747c6ae99SBarry Smith 
94847c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
94947c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
95047c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
95147c6ae99SBarry Smith 
95247c6ae99SBarry Smith   /* determine the matrix preallocation information */
95347c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
95447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
95547c6ae99SBarry Smith 
95647c6ae99SBarry Smith     pstart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
95747c6ae99SBarry Smith     pend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
95847c6ae99SBarry Smith 
95947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
96047c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
96147c6ae99SBarry Smith 
96247c6ae99SBarry Smith       lstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
96347c6ae99SBarry Smith       lend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
96447c6ae99SBarry Smith 
96547c6ae99SBarry Smith       for (k=0; k<nc; k++) {
96647c6ae99SBarry Smith         cnt  = 0;
96747c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
96847c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
96947c6ae99SBarry Smith             if (l || p) {
97047c6ae99SBarry Smith 	      if ((st == DA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
97147c6ae99SBarry Smith                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
97247c6ae99SBarry Smith 		  cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
97347c6ae99SBarry Smith 	      }
97447c6ae99SBarry Smith             } else {
97547c6ae99SBarry Smith 	      if (dfill) {
97647c6ae99SBarry Smith 		for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
97747c6ae99SBarry Smith 		  cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
97847c6ae99SBarry Smith 	      } else {
97947c6ae99SBarry Smith 		for (ifill_col=0; ifill_col<nc; ifill_col++)
98047c6ae99SBarry Smith 		  cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
98147c6ae99SBarry Smith 	      }
98247c6ae99SBarry Smith             }
98347c6ae99SBarry Smith 	  }
98447c6ae99SBarry Smith 	}
98547c6ae99SBarry Smith 	row = k + nc*(slot);
98647c6ae99SBarry Smith         ierr = MatPreallocateSetLocal(ltog,1,&row,cnt,cols,dnz,onz);CHKERRQ(ierr);
98747c6ae99SBarry Smith       }
98847c6ae99SBarry Smith     }
98947c6ae99SBarry Smith   }
99047c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
99147c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
99247c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
99347c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
99447c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
99547c6ae99SBarry Smith 
99647c6ae99SBarry Smith   /*
99747c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
99847c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
99947c6ae99SBarry Smith     PETSc ordering.
100047c6ae99SBarry Smith   */
100147c6ae99SBarry Smith   if (!dd->prealloc_only) {
100247c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
100347c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
100447c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
100547c6ae99SBarry Smith 
100647c6ae99SBarry Smith       pstart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
100747c6ae99SBarry Smith       pend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
100847c6ae99SBarry Smith 
100947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
101047c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
101147c6ae99SBarry Smith 
101247c6ae99SBarry Smith 	lstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
101347c6ae99SBarry Smith 	lend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
101447c6ae99SBarry Smith 
101547c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
101647c6ae99SBarry Smith 	  cnt  = 0;
101747c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
101847c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
101947c6ae99SBarry Smith 	      if (l || p) {
102047c6ae99SBarry Smith 		if ((st == DA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
102147c6ae99SBarry Smith 		  for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
102247c6ae99SBarry Smith 		    cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
102347c6ae99SBarry Smith 		}
102447c6ae99SBarry Smith 	      } else {
102547c6ae99SBarry Smith 		if (dfill) {
102647c6ae99SBarry Smith 		  for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
102747c6ae99SBarry Smith 		    cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
102847c6ae99SBarry Smith 		} else {
102947c6ae99SBarry Smith 		  for (ifill_col=0; ifill_col<nc; ifill_col++)
103047c6ae99SBarry Smith 		    cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
103147c6ae99SBarry Smith 		}
103247c6ae99SBarry Smith 	      }
103347c6ae99SBarry Smith 	    }
103447c6ae99SBarry Smith 	  }
103547c6ae99SBarry Smith 	  row  = k + nc*(slot);
103647c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
103747c6ae99SBarry Smith 	}
103847c6ae99SBarry Smith       }
103947c6ae99SBarry Smith     }
104047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
104147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
104247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
104347c6ae99SBarry Smith   }
104447c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
104547c6ae99SBarry Smith   PetscFunctionReturn(0);
104647c6ae99SBarry Smith }
104747c6ae99SBarry Smith 
104847c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
104947c6ae99SBarry Smith 
105047c6ae99SBarry Smith #undef __FUNCT__
105147c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix3d_MPIAIJ"
1052*9a42bb27SBarry Smith PetscErrorCode DAGetMatrix3d_MPIAIJ(DM da,Mat J)
105347c6ae99SBarry Smith {
105447c6ae99SBarry Smith   PetscErrorCode         ierr;
105547c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
105647c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL;
105747c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
105847c6ae99SBarry Smith   MPI_Comm               comm;
105947c6ae99SBarry Smith   PetscScalar            *values;
106047c6ae99SBarry Smith   DAPeriodicType         wrap;
106147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
106247c6ae99SBarry Smith   DAStencilType          st;
106347c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
106447c6ae99SBarry Smith 
106547c6ae99SBarry Smith   PetscFunctionBegin;
106647c6ae99SBarry Smith   /*
106747c6ae99SBarry Smith          nc - number of components per grid point
106847c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
106947c6ae99SBarry Smith 
107047c6ae99SBarry Smith   */
107147c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
107247c6ae99SBarry Smith   col    = 2*s + 1;
107347c6ae99SBarry Smith 
107447c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
107547c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
107647c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
107747c6ae99SBarry Smith 
107847c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
107947c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
108047c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
108147c6ae99SBarry Smith 
108247c6ae99SBarry Smith   /* determine the matrix preallocation information */
108347c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
108447c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
108547c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
108647c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
108747c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
108847c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
108947c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
109047c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
109147c6ae99SBarry Smith 	kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
109247c6ae99SBarry Smith 	kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
109347c6ae99SBarry Smith 
109447c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
109547c6ae99SBarry Smith 
109647c6ae99SBarry Smith 	cnt  = 0;
109747c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
109847c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
109947c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
110047c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
110147c6ae99SBarry Smith 		if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
110247c6ae99SBarry Smith 		  cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
110347c6ae99SBarry Smith 		}
110447c6ae99SBarry Smith 	      }
110547c6ae99SBarry Smith 	    }
110647c6ae99SBarry Smith 	  }
110747c6ae99SBarry Smith 	  rows[l] = l + nc*(slot);
110847c6ae99SBarry Smith 	}
110947c6ae99SBarry Smith 	ierr = MatPreallocateSetLocal(ltog,nc,rows,cnt,cols,dnz,onz);CHKERRQ(ierr);
111047c6ae99SBarry Smith       }
111147c6ae99SBarry Smith     }
111247c6ae99SBarry Smith   }
111347c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
111447c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
111547c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
111647c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
111747c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
111847c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
111947c6ae99SBarry Smith 
112047c6ae99SBarry Smith   /*
112147c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
112247c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
112347c6ae99SBarry Smith     PETSc ordering.
112447c6ae99SBarry Smith   */
112547c6ae99SBarry Smith   if (!dd->prealloc_only) {
112647c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
112747c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
112847c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
112947c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
113047c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
113147c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
113247c6ae99SBarry Smith 	jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
113347c6ae99SBarry Smith 	jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
113447c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
113547c6ae99SBarry Smith 	  kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
113647c6ae99SBarry Smith 	  kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
113747c6ae99SBarry Smith 
113847c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
113947c6ae99SBarry Smith 
114047c6ae99SBarry Smith 	  cnt  = 0;
114147c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
114247c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
114347c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
114447c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
114547c6ae99SBarry Smith 		  if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
114647c6ae99SBarry Smith 		    cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
114747c6ae99SBarry Smith 		  }
114847c6ae99SBarry Smith 		}
114947c6ae99SBarry Smith 	      }
115047c6ae99SBarry Smith 	    }
115147c6ae99SBarry Smith 	    rows[l]      = l + nc*(slot);
115247c6ae99SBarry Smith 	  }
115347c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
115447c6ae99SBarry Smith 	}
115547c6ae99SBarry Smith       }
115647c6ae99SBarry Smith     }
115747c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
115847c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
115947c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
116047c6ae99SBarry Smith   }
116147c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
116247c6ae99SBarry Smith   PetscFunctionReturn(0);
116347c6ae99SBarry Smith }
116447c6ae99SBarry Smith 
116547c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
116647c6ae99SBarry Smith 
116747c6ae99SBarry Smith #undef __FUNCT__
116847c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix1d_MPIAIJ"
1169*9a42bb27SBarry Smith PetscErrorCode DAGetMatrix1d_MPIAIJ(DM da,Mat J)
117047c6ae99SBarry Smith {
117147c6ae99SBarry Smith   PetscErrorCode         ierr;
117247c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
117347c6ae99SBarry Smith   PetscInt               m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l;
117447c6ae99SBarry Smith   PetscInt               istart,iend;
117547c6ae99SBarry Smith   PetscScalar            *values;
117647c6ae99SBarry Smith   DAPeriodicType         wrap;
117747c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
117847c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
117947c6ae99SBarry Smith 
118047c6ae99SBarry Smith   PetscFunctionBegin;
118147c6ae99SBarry Smith   /*
118247c6ae99SBarry Smith          nc - number of components per grid point
118347c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
118447c6ae99SBarry Smith 
118547c6ae99SBarry Smith   */
118647c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&wrap,0);CHKERRQ(ierr);
118747c6ae99SBarry Smith   col    = 2*s + 1;
118847c6ae99SBarry Smith 
118947c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
119047c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
119147c6ae99SBarry Smith 
119247c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
119347c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
119447c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
119547c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
119647c6ae99SBarry Smith 
119747c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
119847c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
119947c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
120047c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
120147c6ae99SBarry Smith 
120247c6ae99SBarry Smith   /*
120347c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
120447c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
120547c6ae99SBarry Smith     PETSc ordering.
120647c6ae99SBarry Smith   */
120747c6ae99SBarry Smith   if (!dd->prealloc_only) {
120847c6ae99SBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
120947c6ae99SBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
121047c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
121147c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
121247c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
121347c6ae99SBarry Smith       slot   = i - gxs;
121447c6ae99SBarry Smith 
121547c6ae99SBarry Smith       cnt  = 0;
121647c6ae99SBarry Smith       for (l=0; l<nc; l++) {
121747c6ae99SBarry Smith 	for (i1=istart; i1<iend+1; i1++) {
121847c6ae99SBarry Smith 	  cols[cnt++] = l + nc*(slot + i1);
121947c6ae99SBarry Smith 	}
122047c6ae99SBarry Smith 	rows[l]      = l + nc*(slot);
122147c6ae99SBarry Smith       }
122247c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
122347c6ae99SBarry Smith     }
122447c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
122547c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
122647c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
122747c6ae99SBarry Smith   }
122847c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
122947c6ae99SBarry Smith   PetscFunctionReturn(0);
123047c6ae99SBarry Smith }
123147c6ae99SBarry Smith 
123247c6ae99SBarry Smith #undef __FUNCT__
123347c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix2d_MPIBAIJ"
1234*9a42bb27SBarry Smith PetscErrorCode DAGetMatrix2d_MPIBAIJ(DM da,Mat J)
123547c6ae99SBarry Smith {
123647c6ae99SBarry Smith   PetscErrorCode         ierr;
123747c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
123847c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
123947c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
124047c6ae99SBarry Smith   MPI_Comm               comm;
124147c6ae99SBarry Smith   PetscScalar            *values;
124247c6ae99SBarry Smith   DAPeriodicType         wrap;
124347c6ae99SBarry Smith   DAStencilType          st;
124447c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
124547c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
124647c6ae99SBarry Smith 
124747c6ae99SBarry Smith   PetscFunctionBegin;
124847c6ae99SBarry Smith   /*
124947c6ae99SBarry Smith      nc - number of components per grid point
125047c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
125147c6ae99SBarry Smith   */
125247c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
125347c6ae99SBarry Smith   col = 2*s + 1;
125447c6ae99SBarry Smith 
125547c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
125647c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
125747c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
125847c6ae99SBarry Smith 
125947c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
126047c6ae99SBarry Smith 
126147c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
126247c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
126347c6ae99SBarry Smith 
126447c6ae99SBarry Smith   /* determine the matrix preallocation information */
126547c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
126647c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
126747c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
126847c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
126947c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
127047c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
127147c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
127247c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
127347c6ae99SBarry Smith 
127447c6ae99SBarry Smith       /* Find block columns in block row */
127547c6ae99SBarry Smith       cnt  = 0;
127647c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
127747c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
127847c6ae99SBarry Smith           if (st == DA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
127947c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
128047c6ae99SBarry Smith           }
128147c6ae99SBarry Smith         }
128247c6ae99SBarry Smith       }
128347c6ae99SBarry Smith       ierr = MatPreallocateSetLocal(ltogb,1,&slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
128447c6ae99SBarry Smith     }
128547c6ae99SBarry Smith   }
128647c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
128747c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
128847c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
128947c6ae99SBarry Smith 
129047c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
129147c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
129247c6ae99SBarry Smith 
129347c6ae99SBarry Smith   /*
129447c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
129547c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
129647c6ae99SBarry Smith     PETSc ordering.
129747c6ae99SBarry Smith   */
129847c6ae99SBarry Smith   if (!dd->prealloc_only) {
129947c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
130047c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
130147c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
130247c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
130347c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
130447c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
130547c6ae99SBarry Smith 	jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
130647c6ae99SBarry Smith 	jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
130747c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
130847c6ae99SBarry Smith 	cnt  = 0;
130947c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
131047c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
131147c6ae99SBarry Smith             if (st == DA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
131247c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
131347c6ae99SBarry Smith             }
131447c6ae99SBarry Smith           }
131547c6ae99SBarry Smith         }
131647c6ae99SBarry Smith 	ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
131747c6ae99SBarry Smith       }
131847c6ae99SBarry Smith     }
131947c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
132047c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
132147c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
132247c6ae99SBarry Smith   }
132347c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
132447c6ae99SBarry Smith   PetscFunctionReturn(0);
132547c6ae99SBarry Smith }
132647c6ae99SBarry Smith 
132747c6ae99SBarry Smith #undef __FUNCT__
132847c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix3d_MPIBAIJ"
1329*9a42bb27SBarry Smith PetscErrorCode DAGetMatrix3d_MPIBAIJ(DM da,Mat J)
133047c6ae99SBarry Smith {
133147c6ae99SBarry Smith   PetscErrorCode         ierr;
133247c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
133347c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
133447c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
133547c6ae99SBarry Smith   MPI_Comm               comm;
133647c6ae99SBarry Smith   PetscScalar            *values;
133747c6ae99SBarry Smith   DAPeriodicType         wrap;
133847c6ae99SBarry Smith   DAStencilType          st;
133947c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
134047c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
134147c6ae99SBarry Smith 
134247c6ae99SBarry Smith   PetscFunctionBegin;
134347c6ae99SBarry Smith   /*
134447c6ae99SBarry Smith          nc - number of components per grid point
134547c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
134647c6ae99SBarry Smith 
134747c6ae99SBarry Smith   */
134847c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
134947c6ae99SBarry Smith   col    = 2*s + 1;
135047c6ae99SBarry Smith 
135147c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
135247c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
135347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
135447c6ae99SBarry Smith 
135547c6ae99SBarry Smith   ierr  = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
135647c6ae99SBarry Smith 
135747c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
135847c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
135947c6ae99SBarry Smith 
136047c6ae99SBarry Smith   /* determine the matrix preallocation information */
136147c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
136247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
136347c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
136447c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
136547c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
136647c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
136747c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
136847c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
136947c6ae99SBarry Smith 	kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
137047c6ae99SBarry Smith 	kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
137147c6ae99SBarry Smith 
137247c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
137347c6ae99SBarry Smith 
137447c6ae99SBarry Smith 	/* Find block columns in block row */
137547c6ae99SBarry Smith 	cnt  = 0;
137647c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
137747c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
137847c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
137947c6ae99SBarry Smith               if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
138047c6ae99SBarry Smith 		cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
138147c6ae99SBarry Smith 	      }
138247c6ae99SBarry Smith 	    }
138347c6ae99SBarry Smith 	  }
138447c6ae99SBarry Smith 	}
138547c6ae99SBarry Smith 	ierr = MatPreallocateSetLocal(ltogb,1,&slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
138647c6ae99SBarry Smith       }
138747c6ae99SBarry Smith     }
138847c6ae99SBarry Smith   }
138947c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
139047c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
139147c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
139247c6ae99SBarry Smith 
139347c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
139447c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
139547c6ae99SBarry Smith 
139647c6ae99SBarry Smith   /*
139747c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
139847c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
139947c6ae99SBarry Smith     PETSc ordering.
140047c6ae99SBarry Smith   */
140147c6ae99SBarry Smith   if (!dd->prealloc_only) {
140247c6ae99SBarry Smith     ierr  = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
140347c6ae99SBarry Smith     ierr  = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
140447c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
140547c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
140647c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
140747c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
140847c6ae99SBarry Smith 	jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
140947c6ae99SBarry Smith 	jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
141047c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
141147c6ae99SBarry Smith 	  kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
141247c6ae99SBarry Smith 	  kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
141347c6ae99SBarry Smith 
141447c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
141547c6ae99SBarry Smith 
141647c6ae99SBarry Smith 	  cnt  = 0;
141747c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
141847c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
141947c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
142047c6ae99SBarry Smith                 if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
142147c6ae99SBarry Smith                   cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
142247c6ae99SBarry Smith                 }
142347c6ae99SBarry Smith               }
142447c6ae99SBarry Smith             }
142547c6ae99SBarry Smith           }
142647c6ae99SBarry Smith 	  ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
142747c6ae99SBarry Smith 	}
142847c6ae99SBarry Smith       }
142947c6ae99SBarry Smith     }
143047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
143147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
143247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
143347c6ae99SBarry Smith   }
143447c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
143547c6ae99SBarry Smith   PetscFunctionReturn(0);
143647c6ae99SBarry Smith }
143747c6ae99SBarry Smith 
143847c6ae99SBarry Smith #undef __FUNCT__
143947c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular"
144047c6ae99SBarry Smith /*
144147c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
144247c6ae99SBarry Smith   identify in the local ordering with periodic domain.
144347c6ae99SBarry Smith */
144447c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
144547c6ae99SBarry Smith {
144647c6ae99SBarry Smith   PetscErrorCode ierr;
144747c6ae99SBarry Smith   PetscInt i,n;
144847c6ae99SBarry Smith 
144947c6ae99SBarry Smith   PetscFunctionBegin;
145047c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr);
145147c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr);
145247c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
145347c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
145447c6ae99SBarry Smith   }
145547c6ae99SBarry Smith   *cnt = n;
145647c6ae99SBarry Smith   PetscFunctionReturn(0);
145747c6ae99SBarry Smith }
145847c6ae99SBarry Smith 
145947c6ae99SBarry Smith #undef __FUNCT__
146047c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix2d_MPISBAIJ"
1461*9a42bb27SBarry Smith PetscErrorCode DAGetMatrix2d_MPISBAIJ(DM da,Mat J)
146247c6ae99SBarry Smith {
146347c6ae99SBarry Smith   PetscErrorCode         ierr;
146447c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
146547c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
146647c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
146747c6ae99SBarry Smith   MPI_Comm               comm;
146847c6ae99SBarry Smith   PetscScalar            *values;
146947c6ae99SBarry Smith   DAPeriodicType         wrap;
147047c6ae99SBarry Smith   DAStencilType          st;
147147c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
147247c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
147347c6ae99SBarry Smith 
147447c6ae99SBarry Smith   PetscFunctionBegin;
147547c6ae99SBarry Smith   /*
147647c6ae99SBarry Smith      nc - number of components per grid point
147747c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
147847c6ae99SBarry Smith   */
147947c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
148047c6ae99SBarry Smith   col = 2*s + 1;
148147c6ae99SBarry Smith 
148247c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
148347c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
148447c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
148547c6ae99SBarry Smith 
148647c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
148747c6ae99SBarry Smith 
148847c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
148947c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
149047c6ae99SBarry Smith 
149147c6ae99SBarry Smith   /* determine the matrix preallocation information */
149247c6ae99SBarry Smith   ierr = MatPreallocateSymmetricInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
149347c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
149447c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
149547c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
149647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
149747c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
149847c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
149947c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
150047c6ae99SBarry Smith 
150147c6ae99SBarry Smith       /* Find block columns in block row */
150247c6ae99SBarry Smith       cnt  = 0;
150347c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
150447c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
150547c6ae99SBarry Smith           if (st == DA_STENCIL_BOX || !ii || !jj) {
150647c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
150747c6ae99SBarry Smith           }
150847c6ae99SBarry Smith         }
150947c6ae99SBarry Smith       }
151047c6ae99SBarry Smith       ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
151147c6ae99SBarry Smith       ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
151247c6ae99SBarry Smith     }
151347c6ae99SBarry Smith   }
151447c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
151547c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
151647c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
151747c6ae99SBarry Smith 
151847c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
151947c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
152047c6ae99SBarry Smith 
152147c6ae99SBarry Smith   /*
152247c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
152347c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
152447c6ae99SBarry Smith     PETSc ordering.
152547c6ae99SBarry Smith   */
152647c6ae99SBarry Smith   if (!dd->prealloc_only) {
152747c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
152847c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
152947c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
153047c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
153147c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
153247c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
153347c6ae99SBarry Smith         jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
153447c6ae99SBarry Smith         jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
153547c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
153647c6ae99SBarry Smith 
153747c6ae99SBarry Smith         /* Find block columns in block row */
153847c6ae99SBarry Smith         cnt  = 0;
153947c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
154047c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
154147c6ae99SBarry Smith             if (st == DA_STENCIL_BOX || !ii || !jj) {
154247c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
154347c6ae99SBarry Smith             }
154447c6ae99SBarry Smith           }
154547c6ae99SBarry Smith         }
154647c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
154747c6ae99SBarry Smith 	ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
154847c6ae99SBarry Smith       }
154947c6ae99SBarry Smith     }
155047c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
155147c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
155247c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
155347c6ae99SBarry Smith   }
155447c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
155547c6ae99SBarry Smith   PetscFunctionReturn(0);
155647c6ae99SBarry Smith }
155747c6ae99SBarry Smith 
155847c6ae99SBarry Smith #undef __FUNCT__
155947c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix3d_MPISBAIJ"
1560*9a42bb27SBarry Smith PetscErrorCode DAGetMatrix3d_MPISBAIJ(DM da,Mat J)
156147c6ae99SBarry Smith {
156247c6ae99SBarry Smith   PetscErrorCode         ierr;
156347c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
156447c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
156547c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
156647c6ae99SBarry Smith   MPI_Comm               comm;
156747c6ae99SBarry Smith   PetscScalar            *values;
156847c6ae99SBarry Smith   DAPeriodicType         wrap;
156947c6ae99SBarry Smith   DAStencilType          st;
157047c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
157147c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
157247c6ae99SBarry Smith 
157347c6ae99SBarry Smith   PetscFunctionBegin;
157447c6ae99SBarry Smith   /*
157547c6ae99SBarry Smith      nc - number of components per grid point
157647c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
157747c6ae99SBarry Smith   */
157847c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
157947c6ae99SBarry Smith   col = 2*s + 1;
158047c6ae99SBarry Smith 
158147c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
158247c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
158347c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
158447c6ae99SBarry Smith 
158547c6ae99SBarry Smith   /* create the matrix */
158647c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
158747c6ae99SBarry Smith 
158847c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
158947c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
159047c6ae99SBarry Smith 
159147c6ae99SBarry Smith   /* determine the matrix preallocation information */
159247c6ae99SBarry Smith   ierr = MatPreallocateSymmetricInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
159347c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
159447c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
159547c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
159647c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
159747c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
159847c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
159947c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
160047c6ae99SBarry Smith         kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
160147c6ae99SBarry Smith 	kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
160247c6ae99SBarry Smith 
160347c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
160447c6ae99SBarry Smith 
160547c6ae99SBarry Smith 	/* Find block columns in block row */
160647c6ae99SBarry Smith 	cnt  = 0;
160747c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
160847c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
160947c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
161047c6ae99SBarry Smith               if ((st == DA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
161147c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
161247c6ae99SBarry Smith               }
161347c6ae99SBarry Smith             }
161447c6ae99SBarry Smith           }
161547c6ae99SBarry Smith         }
161647c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
161747c6ae99SBarry Smith         ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
161847c6ae99SBarry Smith       }
161947c6ae99SBarry Smith     }
162047c6ae99SBarry Smith   }
162147c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
162247c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
162347c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
162447c6ae99SBarry Smith 
162547c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
162647c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
162747c6ae99SBarry Smith 
162847c6ae99SBarry Smith   /*
162947c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
163047c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
163147c6ae99SBarry Smith     PETSc ordering.
163247c6ae99SBarry Smith   */
163347c6ae99SBarry Smith   if (!dd->prealloc_only) {
163447c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
163547c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
163647c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
163747c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
163847c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
163947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
164047c6ae99SBarry Smith         jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
164147c6ae99SBarry Smith 	jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
164247c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
164347c6ae99SBarry Smith           kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
164447c6ae99SBarry Smith 	  kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
164547c6ae99SBarry Smith 
164647c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
164747c6ae99SBarry Smith 
164847c6ae99SBarry Smith 	  cnt  = 0;
164947c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
165047c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
165147c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
165247c6ae99SBarry Smith                 if ((st == DA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
165347c6ae99SBarry Smith 		  cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
165447c6ae99SBarry Smith 		}
165547c6ae99SBarry Smith 	      }
165647c6ae99SBarry Smith 	    }
165747c6ae99SBarry Smith 	  }
165847c6ae99SBarry Smith           ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
165947c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
166047c6ae99SBarry Smith 	}
166147c6ae99SBarry Smith       }
166247c6ae99SBarry Smith     }
166347c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
166447c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166547c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
166647c6ae99SBarry Smith   }
166747c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
166847c6ae99SBarry Smith   PetscFunctionReturn(0);
166947c6ae99SBarry Smith }
167047c6ae99SBarry Smith 
167147c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
167247c6ae99SBarry Smith 
167347c6ae99SBarry Smith #undef __FUNCT__
167447c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix3d_MPIAIJ_Fill"
1675*9a42bb27SBarry Smith PetscErrorCode DAGetMatrix3d_MPIAIJ_Fill(DM da,Mat J)
167647c6ae99SBarry Smith {
167747c6ae99SBarry Smith   PetscErrorCode         ierr;
167847c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
167947c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz;
168047c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
168147c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
168247c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
168347c6ae99SBarry Smith   MPI_Comm               comm;
168447c6ae99SBarry Smith   PetscScalar            *values;
168547c6ae99SBarry Smith   DAPeriodicType         wrap;
168647c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
168747c6ae99SBarry Smith   DAStencilType          st;
168847c6ae99SBarry Smith 
168947c6ae99SBarry Smith   PetscFunctionBegin;
169047c6ae99SBarry Smith   /*
169147c6ae99SBarry Smith          nc - number of components per grid point
169247c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
169347c6ae99SBarry Smith 
169447c6ae99SBarry Smith   */
169547c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
169647c6ae99SBarry Smith   col    = 2*s + 1;
169747c6ae99SBarry Smith   if (DAXPeriodic(wrap) && (m % col)){
169847c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
169947c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
170047c6ae99SBarry Smith   }
170147c6ae99SBarry Smith   if (DAYPeriodic(wrap) && (n % col)){
170247c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
170347c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
170447c6ae99SBarry Smith   }
170547c6ae99SBarry Smith   if (DAZPeriodic(wrap) && (p % col)){
170647c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
170747c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
170847c6ae99SBarry Smith   }
170947c6ae99SBarry Smith 
171047c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
171147c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
171247c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
171347c6ae99SBarry Smith 
171447c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
171547c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
171647c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
171747c6ae99SBarry Smith 
171847c6ae99SBarry Smith   /* determine the matrix preallocation information */
171947c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
172047c6ae99SBarry Smith 
172147c6ae99SBarry Smith 
172247c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
172347c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
172447c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
172547c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
172647c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
172747c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
172847c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
172947c6ae99SBarry Smith         kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
173047c6ae99SBarry Smith         kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
173147c6ae99SBarry Smith 
173247c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
173347c6ae99SBarry Smith 
173447c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
173547c6ae99SBarry Smith 	  cnt  = 0;
173647c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
173747c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
173847c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
173947c6ae99SBarry Smith 		if (ii || jj || kk) {
174047c6ae99SBarry Smith 		  if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
174147c6ae99SBarry Smith 		    for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
174247c6ae99SBarry Smith 		      cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
174347c6ae99SBarry Smith 		  }
174447c6ae99SBarry Smith 		} else {
174547c6ae99SBarry Smith 		  if (dfill) {
174647c6ae99SBarry Smith 		    for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
174747c6ae99SBarry Smith 		      cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
174847c6ae99SBarry Smith 		  } else {
174947c6ae99SBarry Smith 		    for (ifill_col=0; ifill_col<nc; ifill_col++)
175047c6ae99SBarry Smith 		      cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
175147c6ae99SBarry Smith 		  }
175247c6ae99SBarry Smith 		}
175347c6ae99SBarry Smith 	      }
175447c6ae99SBarry Smith 	    }
175547c6ae99SBarry Smith 	  }
175647c6ae99SBarry Smith 	  row  = l + nc*(slot);
175747c6ae99SBarry Smith 	  ierr = MatPreallocateSetLocal(ltog,1,&row,cnt,cols,dnz,onz);CHKERRQ(ierr);
175847c6ae99SBarry Smith 	}
175947c6ae99SBarry Smith       }
176047c6ae99SBarry Smith     }
176147c6ae99SBarry Smith   }
176247c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
176347c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
176447c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
176547c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
176647c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
176747c6ae99SBarry Smith 
176847c6ae99SBarry Smith   /*
176947c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
177047c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
177147c6ae99SBarry Smith     PETSc ordering.
177247c6ae99SBarry Smith   */
177347c6ae99SBarry Smith   if (!dd->prealloc_only) {
177447c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
177547c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
177647c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
177747c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
177847c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
177947c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
178047c6ae99SBarry Smith 	jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
178147c6ae99SBarry Smith 	jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
178247c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
178347c6ae99SBarry Smith 	  kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
178447c6ae99SBarry Smith 	  kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
178547c6ae99SBarry Smith 
178647c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
178747c6ae99SBarry Smith 
178847c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
178947c6ae99SBarry Smith 	    cnt  = 0;
179047c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
179147c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
179247c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
179347c6ae99SBarry Smith 		  if (ii || jj || kk) {
179447c6ae99SBarry Smith 		    if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
179547c6ae99SBarry Smith 		      for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
179647c6ae99SBarry Smith 			cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
179747c6ae99SBarry Smith 		    }
179847c6ae99SBarry Smith 		  } else {
179947c6ae99SBarry Smith 		    if (dfill) {
180047c6ae99SBarry Smith 		      for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
180147c6ae99SBarry Smith 			cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
180247c6ae99SBarry Smith 		    } else {
180347c6ae99SBarry Smith 		      for (ifill_col=0; ifill_col<nc; ifill_col++)
180447c6ae99SBarry Smith 			cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
180547c6ae99SBarry Smith 		    }
180647c6ae99SBarry Smith 		  }
180747c6ae99SBarry Smith 		}
180847c6ae99SBarry Smith 	      }
180947c6ae99SBarry Smith 	    }
181047c6ae99SBarry Smith 	    row  = l + nc*(slot);
181147c6ae99SBarry Smith 	    ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
181247c6ae99SBarry Smith 	  }
181347c6ae99SBarry Smith 	}
181447c6ae99SBarry Smith       }
181547c6ae99SBarry Smith     }
181647c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
181747c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
181847c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
181947c6ae99SBarry Smith   }
182047c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
182147c6ae99SBarry Smith   PetscFunctionReturn(0);
182247c6ae99SBarry Smith }
1823