xref: /petsc/src/dm/impls/da/fdda.c (revision 47c6ae997ffd1b2afd66b6474dff5950ae8613d1)
1*47c6ae99SBarry Smith #define PETSCDM_DLL
2*47c6ae99SBarry Smith 
3*47c6ae99SBarry Smith #include "private/daimpl.h" /*I      "petscda.h"     I*/
4*47c6ae99SBarry Smith #include "petscmat.h"         /*I      "petscmat.h"    I*/
5*47c6ae99SBarry Smith #include "private/matimpl.h"
6*47c6ae99SBarry Smith 
7*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetColoring1d_MPIAIJ(DA,ISColoringType,ISColoring *);
8*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetColoring2d_MPIAIJ(DA,ISColoringType,ISColoring *);
9*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetColoring2d_5pt_MPIAIJ(DA,ISColoringType,ISColoring *);
10*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetColoring3d_MPIAIJ(DA,ISColoringType,ISColoring *);
11*47c6ae99SBarry Smith 
12*47c6ae99SBarry Smith /*
13*47c6ae99SBarry Smith    For ghost i that may be negative or greater than the upper bound this
14*47c6ae99SBarry Smith   maps it into the 0:m-1 range using periodicity
15*47c6ae99SBarry Smith */
16*47c6ae99SBarry Smith #define SetInRange(i,m) ((i < 0) ? m+i:((i >= m) ? i-m:i))
17*47c6ae99SBarry Smith 
18*47c6ae99SBarry Smith #undef __FUNCT__
19*47c6ae99SBarry Smith #define __FUNCT__ "DASetBlockFills_Private"
20*47c6ae99SBarry Smith static PetscErrorCode DASetBlockFills_Private(PetscInt *dfill,PetscInt w,PetscInt **rfill)
21*47c6ae99SBarry Smith {
22*47c6ae99SBarry Smith   PetscErrorCode ierr;
23*47c6ae99SBarry Smith   PetscInt       i,j,nz,*fill;
24*47c6ae99SBarry Smith 
25*47c6ae99SBarry Smith   PetscFunctionBegin;
26*47c6ae99SBarry Smith   if (!dfill) PetscFunctionReturn(0);
27*47c6ae99SBarry Smith 
28*47c6ae99SBarry Smith   /* count number nonzeros */
29*47c6ae99SBarry Smith   nz = 0;
30*47c6ae99SBarry Smith   for (i=0; i<w; i++) {
31*47c6ae99SBarry Smith     for (j=0; j<w; j++) {
32*47c6ae99SBarry Smith       if (dfill[w*i+j]) nz++;
33*47c6ae99SBarry Smith     }
34*47c6ae99SBarry Smith   }
35*47c6ae99SBarry Smith   ierr = PetscMalloc((nz + w + 1)*sizeof(PetscInt),&fill);CHKERRQ(ierr);
36*47c6ae99SBarry Smith   /* construct modified CSR storage of nonzero structure */
37*47c6ae99SBarry Smith   nz = w + 1;
38*47c6ae99SBarry Smith   for (i=0; i<w; i++) {
39*47c6ae99SBarry Smith     fill[i] = nz;
40*47c6ae99SBarry Smith     for (j=0; j<w; j++) {
41*47c6ae99SBarry Smith       if (dfill[w*i+j]) {
42*47c6ae99SBarry Smith 	fill[nz] = j;
43*47c6ae99SBarry Smith 	nz++;
44*47c6ae99SBarry Smith       }
45*47c6ae99SBarry Smith     }
46*47c6ae99SBarry Smith   }
47*47c6ae99SBarry Smith   fill[w] = nz;
48*47c6ae99SBarry Smith 
49*47c6ae99SBarry Smith   *rfill = fill;
50*47c6ae99SBarry Smith   PetscFunctionReturn(0);
51*47c6ae99SBarry Smith }
52*47c6ae99SBarry Smith 
53*47c6ae99SBarry Smith #undef __FUNCT__
54*47c6ae99SBarry Smith #define __FUNCT__ "DASetMatPreallocateOnly"
55*47c6ae99SBarry Smith /*@
56*47c6ae99SBarry Smith     DASetMatPreallocateOnly - When DAGetMatrix() is called the matrix will be properly
57*47c6ae99SBarry Smith        preallocated but the nonzero structure and zero values will not be set.
58*47c6ae99SBarry Smith 
59*47c6ae99SBarry Smith     Logically Collective on DA
60*47c6ae99SBarry Smith 
61*47c6ae99SBarry Smith     Input Parameter:
62*47c6ae99SBarry Smith +   da - the distributed array
63*47c6ae99SBarry Smith -   only - PETSC_TRUE if only want preallocation
64*47c6ae99SBarry Smith 
65*47c6ae99SBarry Smith 
66*47c6ae99SBarry Smith     Level: developer
67*47c6ae99SBarry Smith 
68*47c6ae99SBarry Smith .seealso DAGetMatrix(), DASetGetMatrix(), DASetBlockSize(), DASetBlockFills()
69*47c6ae99SBarry Smith 
70*47c6ae99SBarry Smith @*/
71*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetMatPreallocateOnly(DA da,PetscBool  only)
72*47c6ae99SBarry Smith {
73*47c6ae99SBarry Smith   DM_DA *dd = (DM_DA*)da->data;
74*47c6ae99SBarry Smith 
75*47c6ae99SBarry Smith   PetscFunctionBegin;
76*47c6ae99SBarry Smith   dd->prealloc_only = only;
77*47c6ae99SBarry Smith   PetscFunctionReturn(0);
78*47c6ae99SBarry Smith }
79*47c6ae99SBarry Smith 
80*47c6ae99SBarry Smith #undef __FUNCT__
81*47c6ae99SBarry Smith #define __FUNCT__ "DASetBlockFills"
82*47c6ae99SBarry Smith /*@
83*47c6ae99SBarry Smith     DASetBlockFills - Sets the fill pattern in each block for a multi-component problem
84*47c6ae99SBarry Smith     of the matrix returned by DAGetMatrix().
85*47c6ae99SBarry Smith 
86*47c6ae99SBarry Smith     Logically Collective on DA
87*47c6ae99SBarry Smith 
88*47c6ae99SBarry Smith     Input Parameter:
89*47c6ae99SBarry Smith +   da - the distributed array
90*47c6ae99SBarry Smith .   dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block)
91*47c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
92*47c6ae99SBarry Smith 
93*47c6ae99SBarry Smith 
94*47c6ae99SBarry Smith     Level: developer
95*47c6ae99SBarry Smith 
96*47c6ae99SBarry Smith     Notes: This only makes sense when you are doing multicomponent problems but using the
97*47c6ae99SBarry Smith        MPIAIJ matrix format
98*47c6ae99SBarry Smith 
99*47c6ae99SBarry Smith            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
100*47c6ae99SBarry Smith        representing coupling and 0 entries for missing coupling. For example
101*47c6ae99SBarry Smith $             dfill[9] = {1, 0, 0,
102*47c6ae99SBarry Smith $                         1, 1, 0,
103*47c6ae99SBarry Smith $                         0, 1, 1}
104*47c6ae99SBarry Smith        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with
105*47c6ae99SBarry Smith        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the
106*47c6ae99SBarry Smith        diagonal block).
107*47c6ae99SBarry Smith 
108*47c6ae99SBarry Smith      DASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
109*47c6ae99SBarry Smith      can be represented in the dfill, ofill format
110*47c6ae99SBarry Smith 
111*47c6ae99SBarry Smith    Contributed by Glenn Hammond
112*47c6ae99SBarry Smith 
113*47c6ae99SBarry Smith .seealso DAGetMatrix(), DASetGetMatrix(), DASetMatPreallocateOnly()
114*47c6ae99SBarry Smith 
115*47c6ae99SBarry Smith @*/
116*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DASetBlockFills(DA da,PetscInt *dfill,PetscInt *ofill)
117*47c6ae99SBarry Smith {
118*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
119*47c6ae99SBarry Smith   PetscErrorCode ierr;
120*47c6ae99SBarry Smith 
121*47c6ae99SBarry Smith   PetscFunctionBegin;
122*47c6ae99SBarry Smith   ierr = DASetBlockFills_Private(dfill,dd->w,&dd->dfill);CHKERRQ(ierr);
123*47c6ae99SBarry Smith   ierr = DASetBlockFills_Private(ofill,dd->w,&dd->ofill);CHKERRQ(ierr);
124*47c6ae99SBarry Smith   PetscFunctionReturn(0);
125*47c6ae99SBarry Smith }
126*47c6ae99SBarry Smith 
127*47c6ae99SBarry Smith 
128*47c6ae99SBarry Smith #undef __FUNCT__
129*47c6ae99SBarry Smith #define __FUNCT__ "DAGetColoring"
130*47c6ae99SBarry Smith /*@
131*47c6ae99SBarry Smith     DAGetColoring - Gets the coloring required for computing the Jacobian via
132*47c6ae99SBarry Smith     finite differences on a function defined using a stencil on the DA.
133*47c6ae99SBarry Smith 
134*47c6ae99SBarry Smith     Collective on DA
135*47c6ae99SBarry Smith 
136*47c6ae99SBarry Smith     Input Parameter:
137*47c6ae99SBarry Smith +   da - the distributed array
138*47c6ae99SBarry Smith .   ctype - IS_COLORING_GLOBAL or IS_COLORING_GHOSTED
139*47c6ae99SBarry Smith -   mtype - either MATAIJ or MATBAIJ
140*47c6ae99SBarry Smith 
141*47c6ae99SBarry Smith     Output Parameters:
142*47c6ae99SBarry Smith .   coloring - matrix coloring for use in computing Jacobians (or PETSC_NULL if not needed)
143*47c6ae99SBarry Smith 
144*47c6ae99SBarry Smith     Level: advanced
145*47c6ae99SBarry Smith 
146*47c6ae99SBarry Smith     Notes: These compute the graph coloring of the graph of A^{T}A. The coloring used
147*47c6ae99SBarry Smith    for efficient (parallel or thread based) triangular solves etc is NOT
148*47c6ae99SBarry Smith    available.
149*47c6ae99SBarry Smith 
150*47c6ae99SBarry Smith         For BAIJ matrices this colors the graph for the blocks, not for the individual matrix elements;
151*47c6ae99SBarry Smith     the same as MatGetColoring().
152*47c6ae99SBarry Smith 
153*47c6ae99SBarry Smith .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), ISColoringType, ISColoring
154*47c6ae99SBarry Smith 
155*47c6ae99SBarry Smith @*/
156*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetColoring(DA da,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
157*47c6ae99SBarry Smith {
158*47c6ae99SBarry Smith   PetscErrorCode ierr;
159*47c6ae99SBarry Smith   PetscInt       dim,m,n,p,nc;
160*47c6ae99SBarry Smith   DAPeriodicType wrap;
161*47c6ae99SBarry Smith   MPI_Comm       comm;
162*47c6ae99SBarry Smith   PetscMPIInt    size;
163*47c6ae99SBarry Smith   PetscBool      isBAIJ;
164*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
165*47c6ae99SBarry Smith 
166*47c6ae99SBarry Smith   PetscFunctionBegin;
167*47c6ae99SBarry Smith   /*
168*47c6ae99SBarry Smith                                   m
169*47c6ae99SBarry Smith           ------------------------------------------------------
170*47c6ae99SBarry Smith          |                                                     |
171*47c6ae99SBarry Smith          |                                                     |
172*47c6ae99SBarry Smith          |               ----------------------                |
173*47c6ae99SBarry Smith          |               |                    |                |
174*47c6ae99SBarry Smith       n  |           yn  |                    |                |
175*47c6ae99SBarry Smith          |               |                    |                |
176*47c6ae99SBarry Smith          |               .---------------------                |
177*47c6ae99SBarry Smith          |             (xs,ys)     xn                          |
178*47c6ae99SBarry Smith          |            .                                        |
179*47c6ae99SBarry Smith          |         (gxs,gys)                                   |
180*47c6ae99SBarry Smith          |                                                     |
181*47c6ae99SBarry Smith           -----------------------------------------------------
182*47c6ae99SBarry Smith   */
183*47c6ae99SBarry Smith 
184*47c6ae99SBarry Smith   /*
185*47c6ae99SBarry Smith          nc - number of components per grid point
186*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
187*47c6ae99SBarry Smith 
188*47c6ae99SBarry Smith   */
189*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&wrap,0);CHKERRQ(ierr);
190*47c6ae99SBarry Smith 
191*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
192*47c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
193*47c6ae99SBarry Smith   if (ctype == IS_COLORING_GHOSTED){
194*47c6ae99SBarry Smith     if (size == 1) {
195*47c6ae99SBarry Smith       ctype = IS_COLORING_GLOBAL;
196*47c6ae99SBarry Smith     } else if (dim > 1){
197*47c6ae99SBarry Smith       if ((m==1 && DAXPeriodic(wrap)) || (n==1 && DAYPeriodic(wrap)) || (p==1 && DAZPeriodic(wrap))){
198*47c6ae99SBarry 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");
199*47c6ae99SBarry Smith       }
200*47c6ae99SBarry Smith     }
201*47c6ae99SBarry Smith   }
202*47c6ae99SBarry Smith 
203*47c6ae99SBarry Smith   /* Tell the DA it has 1 degree of freedom per grid point so that the coloring for BAIJ
204*47c6ae99SBarry Smith      matrices is for the blocks, not the individual matrix elements  */
205*47c6ae99SBarry Smith   ierr = PetscStrcmp(mtype,MATBAIJ,&isBAIJ);CHKERRQ(ierr);
206*47c6ae99SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATMPIBAIJ,&isBAIJ);CHKERRQ(ierr);}
207*47c6ae99SBarry Smith   if (!isBAIJ) {ierr = PetscStrcmp(mtype,MATSEQBAIJ,&isBAIJ);CHKERRQ(ierr);}
208*47c6ae99SBarry Smith   if (isBAIJ) {
209*47c6ae99SBarry Smith     dd->w = 1;
210*47c6ae99SBarry Smith     dd->xs = dd->xs/nc;
211*47c6ae99SBarry Smith     dd->xe = dd->xe/nc;
212*47c6ae99SBarry Smith     dd->Xs = dd->Xs/nc;
213*47c6ae99SBarry Smith     dd->Xe = dd->Xe/nc;
214*47c6ae99SBarry Smith   }
215*47c6ae99SBarry Smith 
216*47c6ae99SBarry Smith   /*
217*47c6ae99SBarry Smith      We do not provide a getcoloring function in the DA operations because
218*47c6ae99SBarry Smith    the basic DA does not know about matrices. We think of DA as being more
219*47c6ae99SBarry Smith    more low-level then matrices.
220*47c6ae99SBarry Smith   */
221*47c6ae99SBarry Smith   if (dim == 1) {
222*47c6ae99SBarry Smith     ierr = DAGetColoring1d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
223*47c6ae99SBarry Smith   } else if (dim == 2) {
224*47c6ae99SBarry Smith     ierr =  DAGetColoring2d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
225*47c6ae99SBarry Smith   } else if (dim == 3) {
226*47c6ae99SBarry Smith     ierr =  DAGetColoring3d_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
227*47c6ae99SBarry Smith   } else {
228*47c6ae99SBarry Smith     SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim);
229*47c6ae99SBarry Smith   }
230*47c6ae99SBarry Smith   if (isBAIJ) {
231*47c6ae99SBarry Smith     dd->w = nc;
232*47c6ae99SBarry Smith     dd->xs = dd->xs*nc;
233*47c6ae99SBarry Smith     dd->xe = dd->xe*nc;
234*47c6ae99SBarry Smith     dd->Xs = dd->Xs*nc;
235*47c6ae99SBarry Smith     dd->Xe = dd->Xe*nc;
236*47c6ae99SBarry Smith   }
237*47c6ae99SBarry Smith   PetscFunctionReturn(0);
238*47c6ae99SBarry Smith }
239*47c6ae99SBarry Smith 
240*47c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
241*47c6ae99SBarry Smith 
242*47c6ae99SBarry Smith #undef __FUNCT__
243*47c6ae99SBarry Smith #define __FUNCT__ "DAGetColoring2d_MPIAIJ"
244*47c6ae99SBarry Smith PetscErrorCode DAGetColoring2d_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring)
245*47c6ae99SBarry Smith {
246*47c6ae99SBarry Smith   PetscErrorCode         ierr;
247*47c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col;
248*47c6ae99SBarry Smith   PetscInt               ncolors;
249*47c6ae99SBarry Smith   MPI_Comm               comm;
250*47c6ae99SBarry Smith   DAPeriodicType         wrap;
251*47c6ae99SBarry Smith   DAStencilType          st;
252*47c6ae99SBarry Smith   ISColoringValue        *colors;
253*47c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
254*47c6ae99SBarry Smith 
255*47c6ae99SBarry Smith   PetscFunctionBegin;
256*47c6ae99SBarry Smith   /*
257*47c6ae99SBarry Smith          nc - number of components per grid point
258*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
259*47c6ae99SBarry Smith 
260*47c6ae99SBarry Smith   */
261*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
262*47c6ae99SBarry Smith   col    = 2*s + 1;
263*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
264*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
265*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
266*47c6ae99SBarry Smith 
267*47c6ae99SBarry Smith   /* special case as taught to us by Paul Hovland */
268*47c6ae99SBarry Smith   if (st == DA_STENCIL_STAR && s == 1) {
269*47c6ae99SBarry Smith     ierr = DAGetColoring2d_5pt_MPIAIJ(da,ctype,coloring);CHKERRQ(ierr);
270*47c6ae99SBarry Smith   } else {
271*47c6ae99SBarry Smith 
272*47c6ae99SBarry Smith     if (DAXPeriodic(wrap) && (m % col)){
273*47c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\
274*47c6ae99SBarry Smith                  by 2*stencil_width + 1 (%d)\n", m, col);
275*47c6ae99SBarry Smith     }
276*47c6ae99SBarry Smith     if (DAYPeriodic(wrap) && (n % col)){
277*47c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\
278*47c6ae99SBarry Smith                  by 2*stencil_width + 1 (%d)\n", n, col);
279*47c6ae99SBarry Smith     }
280*47c6ae99SBarry Smith     if (ctype == IS_COLORING_GLOBAL) {
281*47c6ae99SBarry Smith       if (!dd->localcoloring) {
282*47c6ae99SBarry Smith 	ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
283*47c6ae99SBarry Smith 	ii = 0;
284*47c6ae99SBarry Smith 	for (j=ys; j<ys+ny; j++) {
285*47c6ae99SBarry Smith 	  for (i=xs; i<xs+nx; i++) {
286*47c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
287*47c6ae99SBarry Smith 	      colors[ii++] = k + nc*((i % col) + col*(j % col));
288*47c6ae99SBarry Smith 	    }
289*47c6ae99SBarry Smith 	  }
290*47c6ae99SBarry Smith 	}
291*47c6ae99SBarry Smith         ncolors = nc + nc*(col-1 + col*(col-1));
292*47c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
293*47c6ae99SBarry Smith       }
294*47c6ae99SBarry Smith       *coloring = dd->localcoloring;
295*47c6ae99SBarry Smith     } else if (ctype == IS_COLORING_GHOSTED) {
296*47c6ae99SBarry Smith       if (!dd->ghostedcoloring) {
297*47c6ae99SBarry Smith 	ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
298*47c6ae99SBarry Smith 	ii = 0;
299*47c6ae99SBarry Smith 	for (j=gys; j<gys+gny; j++) {
300*47c6ae99SBarry Smith 	  for (i=gxs; i<gxs+gnx; i++) {
301*47c6ae99SBarry Smith 	    for (k=0; k<nc; k++) {
302*47c6ae99SBarry Smith 	      /* the complicated stuff is to handle periodic boundaries */
303*47c6ae99SBarry Smith 	      colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
304*47c6ae99SBarry Smith 	    }
305*47c6ae99SBarry Smith 	  }
306*47c6ae99SBarry Smith 	}
307*47c6ae99SBarry Smith         ncolors = nc + nc*(col - 1 + col*(col-1));
308*47c6ae99SBarry Smith 	ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
309*47c6ae99SBarry Smith         /* PetscIntView(ncolors,(PetscInt *)colors,0); */
310*47c6ae99SBarry Smith 
311*47c6ae99SBarry Smith 	ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
312*47c6ae99SBarry Smith       }
313*47c6ae99SBarry Smith       *coloring = dd->ghostedcoloring;
314*47c6ae99SBarry Smith     } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
315*47c6ae99SBarry Smith   }
316*47c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
317*47c6ae99SBarry Smith   PetscFunctionReturn(0);
318*47c6ae99SBarry Smith }
319*47c6ae99SBarry Smith 
320*47c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
321*47c6ae99SBarry Smith 
322*47c6ae99SBarry Smith #undef __FUNCT__
323*47c6ae99SBarry Smith #define __FUNCT__ "DAGetColoring3d_MPIAIJ"
324*47c6ae99SBarry Smith PetscErrorCode DAGetColoring3d_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring)
325*47c6ae99SBarry Smith {
326*47c6ae99SBarry Smith   PetscErrorCode  ierr;
327*47c6ae99SBarry 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;
328*47c6ae99SBarry Smith   PetscInt        ncolors;
329*47c6ae99SBarry Smith   MPI_Comm        comm;
330*47c6ae99SBarry Smith   DAPeriodicType  wrap;
331*47c6ae99SBarry Smith   DAStencilType   st;
332*47c6ae99SBarry Smith   ISColoringValue *colors;
333*47c6ae99SBarry Smith   DM_DA           *dd = (DM_DA*)da->data;
334*47c6ae99SBarry Smith 
335*47c6ae99SBarry Smith   PetscFunctionBegin;
336*47c6ae99SBarry Smith   /*
337*47c6ae99SBarry Smith          nc - number of components per grid point
338*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
339*47c6ae99SBarry Smith 
340*47c6ae99SBarry Smith   */
341*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&wrap,&st);CHKERRQ(ierr);
342*47c6ae99SBarry Smith   col    = 2*s + 1;
343*47c6ae99SBarry Smith   if (DAXPeriodic(wrap) && (m % col)){
344*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
345*47c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
346*47c6ae99SBarry Smith   }
347*47c6ae99SBarry Smith   if (DAYPeriodic(wrap) && (n % col)){
348*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
349*47c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
350*47c6ae99SBarry Smith   }
351*47c6ae99SBarry Smith   if (DAZPeriodic(wrap) && (p % col)){
352*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
353*47c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
354*47c6ae99SBarry Smith   }
355*47c6ae99SBarry Smith 
356*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
357*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
358*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
359*47c6ae99SBarry Smith 
360*47c6ae99SBarry Smith   /* create the coloring */
361*47c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
362*47c6ae99SBarry Smith     if (!dd->localcoloring) {
363*47c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
364*47c6ae99SBarry Smith       ii = 0;
365*47c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
366*47c6ae99SBarry Smith         for (j=ys; j<ys+ny; j++) {
367*47c6ae99SBarry Smith           for (i=xs; i<xs+nx; i++) {
368*47c6ae99SBarry Smith             for (l=0; l<nc; l++) {
369*47c6ae99SBarry Smith               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
370*47c6ae99SBarry Smith             }
371*47c6ae99SBarry Smith           }
372*47c6ae99SBarry Smith         }
373*47c6ae99SBarry Smith       }
374*47c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
375*47c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);CHKERRQ(ierr);
376*47c6ae99SBarry Smith     }
377*47c6ae99SBarry Smith     *coloring = dd->localcoloring;
378*47c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
379*47c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
380*47c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
381*47c6ae99SBarry Smith       ii = 0;
382*47c6ae99SBarry Smith       for (k=gzs; k<gzs+gnz; k++) {
383*47c6ae99SBarry Smith         for (j=gys; j<gys+gny; j++) {
384*47c6ae99SBarry Smith           for (i=gxs; i<gxs+gnx; i++) {
385*47c6ae99SBarry Smith             for (l=0; l<nc; l++) {
386*47c6ae99SBarry Smith               /* the complicated stuff is to handle periodic boundaries */
387*47c6ae99SBarry Smith               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
388*47c6ae99SBarry Smith             }
389*47c6ae99SBarry Smith           }
390*47c6ae99SBarry Smith         }
391*47c6ae99SBarry Smith       }
392*47c6ae99SBarry Smith       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
393*47c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
394*47c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
395*47c6ae99SBarry Smith     }
396*47c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
397*47c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
398*47c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
399*47c6ae99SBarry Smith   PetscFunctionReturn(0);
400*47c6ae99SBarry Smith }
401*47c6ae99SBarry Smith 
402*47c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
403*47c6ae99SBarry Smith 
404*47c6ae99SBarry Smith #undef __FUNCT__
405*47c6ae99SBarry Smith #define __FUNCT__ "DAGetColoring1d_MPIAIJ"
406*47c6ae99SBarry Smith PetscErrorCode DAGetColoring1d_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring)
407*47c6ae99SBarry Smith {
408*47c6ae99SBarry Smith   PetscErrorCode  ierr;
409*47c6ae99SBarry Smith   PetscInt        xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
410*47c6ae99SBarry Smith   PetscInt        ncolors;
411*47c6ae99SBarry Smith   MPI_Comm        comm;
412*47c6ae99SBarry Smith   DAPeriodicType  wrap;
413*47c6ae99SBarry Smith   ISColoringValue *colors;
414*47c6ae99SBarry Smith   DM_DA           *dd = (DM_DA*)da->data;
415*47c6ae99SBarry Smith 
416*47c6ae99SBarry Smith   PetscFunctionBegin;
417*47c6ae99SBarry Smith   /*
418*47c6ae99SBarry Smith          nc - number of components per grid point
419*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
420*47c6ae99SBarry Smith 
421*47c6ae99SBarry Smith   */
422*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&wrap,0);CHKERRQ(ierr);
423*47c6ae99SBarry Smith   col    = 2*s + 1;
424*47c6ae99SBarry Smith 
425*47c6ae99SBarry Smith   if (DAXPeriodic(wrap) && (m % col)) {
426*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points is divisible\n\
427*47c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
428*47c6ae99SBarry Smith   }
429*47c6ae99SBarry Smith 
430*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
431*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
432*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
433*47c6ae99SBarry Smith 
434*47c6ae99SBarry Smith   /* create the coloring */
435*47c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
436*47c6ae99SBarry Smith     if (!dd->localcoloring) {
437*47c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
438*47c6ae99SBarry Smith       i1 = 0;
439*47c6ae99SBarry Smith       for (i=xs; i<xs+nx; i++) {
440*47c6ae99SBarry Smith         for (l=0; l<nc; l++) {
441*47c6ae99SBarry Smith           colors[i1++] = l + nc*(i % col);
442*47c6ae99SBarry Smith         }
443*47c6ae99SBarry Smith       }
444*47c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
445*47c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);CHKERRQ(ierr);
446*47c6ae99SBarry Smith     }
447*47c6ae99SBarry Smith     *coloring = dd->localcoloring;
448*47c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
449*47c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
450*47c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
451*47c6ae99SBarry Smith       i1 = 0;
452*47c6ae99SBarry Smith       for (i=gxs; i<gxs+gnx; i++) {
453*47c6ae99SBarry Smith         for (l=0; l<nc; l++) {
454*47c6ae99SBarry Smith           /* the complicated stuff is to handle periodic boundaries */
455*47c6ae99SBarry Smith           colors[i1++] = l + nc*(SetInRange(i,m) % col);
456*47c6ae99SBarry Smith         }
457*47c6ae99SBarry Smith       }
458*47c6ae99SBarry Smith       ncolors = nc + nc*(col-1);
459*47c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
460*47c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
461*47c6ae99SBarry Smith     }
462*47c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
463*47c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
464*47c6ae99SBarry Smith   ierr = ISColoringReference(*coloring);CHKERRQ(ierr);
465*47c6ae99SBarry Smith   PetscFunctionReturn(0);
466*47c6ae99SBarry Smith }
467*47c6ae99SBarry Smith 
468*47c6ae99SBarry Smith #undef __FUNCT__
469*47c6ae99SBarry Smith #define __FUNCT__ "DAGetColoring2d_5pt_MPIAIJ"
470*47c6ae99SBarry Smith PetscErrorCode DAGetColoring2d_5pt_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring)
471*47c6ae99SBarry Smith {
472*47c6ae99SBarry Smith   PetscErrorCode  ierr;
473*47c6ae99SBarry Smith   PetscInt        xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
474*47c6ae99SBarry Smith   PetscInt        ncolors;
475*47c6ae99SBarry Smith   MPI_Comm        comm;
476*47c6ae99SBarry Smith   DAPeriodicType  wrap;
477*47c6ae99SBarry Smith   ISColoringValue *colors;
478*47c6ae99SBarry Smith   DM_DA           *dd = (DM_DA*)da->data;
479*47c6ae99SBarry Smith 
480*47c6ae99SBarry Smith   PetscFunctionBegin;
481*47c6ae99SBarry Smith   /*
482*47c6ae99SBarry Smith          nc - number of components per grid point
483*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
484*47c6ae99SBarry Smith 
485*47c6ae99SBarry Smith   */
486*47c6ae99SBarry Smith   ierr   = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,0);CHKERRQ(ierr);
487*47c6ae99SBarry Smith   ierr   = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
488*47c6ae99SBarry Smith   ierr   = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
489*47c6ae99SBarry Smith   ierr   = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
490*47c6ae99SBarry Smith 
491*47c6ae99SBarry Smith   if (DAXPeriodic(wrap) && (m % 5)){
492*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
493*47c6ae99SBarry Smith                  by 5\n");
494*47c6ae99SBarry Smith   }
495*47c6ae99SBarry Smith   if (DAYPeriodic(wrap) && (n % 5)){
496*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
497*47c6ae99SBarry Smith                  by 5\n");
498*47c6ae99SBarry Smith   }
499*47c6ae99SBarry Smith 
500*47c6ae99SBarry Smith   /* create the coloring */
501*47c6ae99SBarry Smith   if (ctype == IS_COLORING_GLOBAL) {
502*47c6ae99SBarry Smith     if (!dd->localcoloring) {
503*47c6ae99SBarry Smith       ierr = PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
504*47c6ae99SBarry Smith       ii = 0;
505*47c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
506*47c6ae99SBarry Smith 	for (i=xs; i<xs+nx; i++) {
507*47c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
508*47c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*j+i) % 5);
509*47c6ae99SBarry Smith 	  }
510*47c6ae99SBarry Smith 	}
511*47c6ae99SBarry Smith       }
512*47c6ae99SBarry Smith       ncolors = 5*nc;
513*47c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);CHKERRQ(ierr);
514*47c6ae99SBarry Smith     }
515*47c6ae99SBarry Smith     *coloring = dd->localcoloring;
516*47c6ae99SBarry Smith   } else if (ctype == IS_COLORING_GHOSTED) {
517*47c6ae99SBarry Smith     if (!dd->ghostedcoloring) {
518*47c6ae99SBarry Smith       ierr = PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);CHKERRQ(ierr);
519*47c6ae99SBarry Smith       ii = 0;
520*47c6ae99SBarry Smith       for (j=gys; j<gys+gny; j++) {
521*47c6ae99SBarry Smith 	for (i=gxs; i<gxs+gnx; i++) {
522*47c6ae99SBarry Smith 	  for (k=0; k<nc; k++) {
523*47c6ae99SBarry Smith 	    colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
524*47c6ae99SBarry Smith 	  }
525*47c6ae99SBarry Smith 	}
526*47c6ae99SBarry Smith       }
527*47c6ae99SBarry Smith       ncolors = 5*nc;
528*47c6ae99SBarry Smith       ierr = ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);CHKERRQ(ierr);
529*47c6ae99SBarry Smith       ierr = ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);CHKERRQ(ierr);
530*47c6ae99SBarry Smith     }
531*47c6ae99SBarry Smith     *coloring = dd->ghostedcoloring;
532*47c6ae99SBarry Smith   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
533*47c6ae99SBarry Smith   PetscFunctionReturn(0);
534*47c6ae99SBarry Smith }
535*47c6ae99SBarry Smith 
536*47c6ae99SBarry Smith /* =========================================================================== */
537*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetMatrix1d_MPIAIJ(DA,Mat);
538*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetMatrix2d_MPIAIJ(DA,Mat);
539*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetMatrix2d_MPIAIJ_Fill(DA,Mat);
540*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetMatrix3d_MPIAIJ(DA,Mat);
541*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetMatrix3d_MPIAIJ_Fill(DA,Mat);
542*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetMatrix2d_MPIBAIJ(DA,Mat);
543*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetMatrix3d_MPIBAIJ(DA,Mat);
544*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetMatrix2d_MPISBAIJ(DA,Mat);
545*47c6ae99SBarry Smith EXTERN PetscErrorCode DAGetMatrix3d_MPISBAIJ(DA,Mat);
546*47c6ae99SBarry Smith 
547*47c6ae99SBarry Smith #undef __FUNCT__
548*47c6ae99SBarry Smith #define __FUNCT__ "MatSetDA"
549*47c6ae99SBarry Smith /*@
550*47c6ae99SBarry Smith    MatSetDA - Sets the DA that is to be used by the HYPRE_StructMatrix PETSc matrix
551*47c6ae99SBarry Smith 
552*47c6ae99SBarry Smith    Logically Collective on Mat
553*47c6ae99SBarry Smith 
554*47c6ae99SBarry Smith    Input Parameters:
555*47c6ae99SBarry Smith +  mat - the matrix
556*47c6ae99SBarry Smith -  da - the da
557*47c6ae99SBarry Smith 
558*47c6ae99SBarry Smith    Level: intermediate
559*47c6ae99SBarry Smith 
560*47c6ae99SBarry Smith @*/
561*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT MatSetDA(Mat mat,DA da)
562*47c6ae99SBarry Smith {
563*47c6ae99SBarry Smith   PetscErrorCode ierr;
564*47c6ae99SBarry Smith 
565*47c6ae99SBarry Smith   PetscFunctionBegin;
566*47c6ae99SBarry Smith   PetscValidHeaderSpecific(mat,MAT_CLASSID,1);
567*47c6ae99SBarry Smith   PetscValidHeaderSpecific(da,DM_CLASSID,1);
568*47c6ae99SBarry Smith   ierr = PetscTryMethod(mat,"MatSetDA_C",(Mat,DA),(mat,da));CHKERRQ(ierr);
569*47c6ae99SBarry Smith   PetscFunctionReturn(0);
570*47c6ae99SBarry Smith }
571*47c6ae99SBarry Smith 
572*47c6ae99SBarry Smith EXTERN_C_BEGIN
573*47c6ae99SBarry Smith #undef __FUNCT__
574*47c6ae99SBarry Smith #define __FUNCT__ "MatView_MPI_DA"
575*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT MatView_MPI_DA(Mat A,PetscViewer viewer)
576*47c6ae99SBarry Smith {
577*47c6ae99SBarry Smith   DA             da;
578*47c6ae99SBarry Smith   PetscErrorCode ierr;
579*47c6ae99SBarry Smith   const char     *prefix;
580*47c6ae99SBarry Smith   Mat            Anatural;
581*47c6ae99SBarry Smith   AO             ao;
582*47c6ae99SBarry Smith   PetscInt       rstart,rend,*petsc,i;
583*47c6ae99SBarry Smith   IS             is;
584*47c6ae99SBarry Smith   MPI_Comm       comm;
585*47c6ae99SBarry Smith 
586*47c6ae99SBarry Smith   PetscFunctionBegin;
587*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
588*47c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"DA",(PetscObject*)&da);CHKERRQ(ierr);
589*47c6ae99SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DA");
590*47c6ae99SBarry Smith 
591*47c6ae99SBarry Smith   ierr = DAGetAO(da,&ao);CHKERRQ(ierr);
592*47c6ae99SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
593*47c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);CHKERRQ(ierr);
594*47c6ae99SBarry Smith   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
595*47c6ae99SBarry Smith   ierr = AOApplicationToPetsc(ao,rend-rstart,petsc);CHKERRQ(ierr);
596*47c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
597*47c6ae99SBarry Smith 
598*47c6ae99SBarry Smith   /* call viewer on natural ordering */
599*47c6ae99SBarry Smith   ierr = MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);CHKERRQ(ierr);
600*47c6ae99SBarry Smith   ierr = ISDestroy(is);CHKERRQ(ierr);
601*47c6ae99SBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);CHKERRQ(ierr);
602*47c6ae99SBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);CHKERRQ(ierr);
603*47c6ae99SBarry Smith   ierr = PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);CHKERRQ(ierr);
604*47c6ae99SBarry Smith   ierr = MatView(Anatural,viewer);CHKERRQ(ierr);
605*47c6ae99SBarry Smith   ierr = MatDestroy(Anatural);CHKERRQ(ierr);
606*47c6ae99SBarry Smith   PetscFunctionReturn(0);
607*47c6ae99SBarry Smith }
608*47c6ae99SBarry Smith EXTERN_C_END
609*47c6ae99SBarry Smith 
610*47c6ae99SBarry Smith EXTERN_C_BEGIN
611*47c6ae99SBarry Smith #undef __FUNCT__
612*47c6ae99SBarry Smith #define __FUNCT__ "MatLoad_MPI_DA"
613*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT MatLoad_MPI_DA(Mat A,PetscViewer viewer)
614*47c6ae99SBarry Smith {
615*47c6ae99SBarry Smith   DA             da;
616*47c6ae99SBarry Smith   PetscErrorCode ierr;
617*47c6ae99SBarry Smith   Mat            Anatural,Aapp;
618*47c6ae99SBarry Smith   AO             ao;
619*47c6ae99SBarry Smith   PetscInt       rstart,rend,*app,i;
620*47c6ae99SBarry Smith   IS             is;
621*47c6ae99SBarry Smith   MPI_Comm       comm;
622*47c6ae99SBarry Smith 
623*47c6ae99SBarry Smith   PetscFunctionBegin;
624*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
625*47c6ae99SBarry Smith   ierr = PetscObjectQuery((PetscObject)A,"DA",(PetscObject*)&da);CHKERRQ(ierr);
626*47c6ae99SBarry Smith   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DA");
627*47c6ae99SBarry Smith 
628*47c6ae99SBarry Smith   /* Load the matrix in natural ordering */
629*47c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)A)->comm,&Anatural);CHKERRQ(ierr);
630*47c6ae99SBarry Smith   ierr = MatSetType(Anatural,((PetscObject)A)->type_name);CHKERRQ(ierr);
631*47c6ae99SBarry Smith   ierr = MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);CHKERRQ(ierr);
632*47c6ae99SBarry Smith   ierr = MatLoad(Anatural,viewer);CHKERRQ(ierr);
633*47c6ae99SBarry Smith 
634*47c6ae99SBarry Smith   /* Map natural ordering to application ordering and create IS */
635*47c6ae99SBarry Smith   ierr = DAGetAO(da,&ao);CHKERRQ(ierr);
636*47c6ae99SBarry Smith   ierr = MatGetOwnershipRange(Anatural,&rstart,&rend);CHKERRQ(ierr);
637*47c6ae99SBarry Smith   ierr = PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);CHKERRQ(ierr);
638*47c6ae99SBarry Smith   for (i=rstart; i<rend; i++) app[i-rstart] = i;
639*47c6ae99SBarry Smith   ierr = AOPetscToApplication(ao,rend-rstart,app);CHKERRQ(ierr);
640*47c6ae99SBarry Smith   ierr = ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);CHKERRQ(ierr);
641*47c6ae99SBarry Smith 
642*47c6ae99SBarry Smith   /* Do permutation and replace header */
643*47c6ae99SBarry Smith   ierr = MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);CHKERRQ(ierr);
644*47c6ae99SBarry Smith   ierr = MatHeaderReplace(A,Aapp);CHKERRQ(ierr);
645*47c6ae99SBarry Smith   ierr = ISDestroy(is);CHKERRQ(ierr);
646*47c6ae99SBarry Smith   ierr = MatDestroy(Anatural);CHKERRQ(ierr);
647*47c6ae99SBarry Smith   PetscFunctionReturn(0);
648*47c6ae99SBarry Smith }
649*47c6ae99SBarry Smith EXTERN_C_END
650*47c6ae99SBarry Smith 
651*47c6ae99SBarry Smith 
652*47c6ae99SBarry Smith #undef __FUNCT__
653*47c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix"
654*47c6ae99SBarry Smith /*@C
655*47c6ae99SBarry Smith     DAGetMatrix - Creates a matrix with the correct parallel layout and nonzero structure required for
656*47c6ae99SBarry Smith       computing the Jacobian on a function defined using the stencil set in the DA.
657*47c6ae99SBarry Smith 
658*47c6ae99SBarry Smith     Collective on DA
659*47c6ae99SBarry Smith 
660*47c6ae99SBarry Smith     Input Parameter:
661*47c6ae99SBarry Smith +   da - the distributed array
662*47c6ae99SBarry Smith -   mtype - Supported types are MATSEQAIJ, MATMPIAIJ, MATSEQBAIJ, MATMPIBAIJ, MATSEQSBAIJ, MATMPISBAIJ,
663*47c6ae99SBarry Smith             or any type which inherits from one of these (such as MATAIJ, MATBAIJ, MATSBAIJ)
664*47c6ae99SBarry Smith 
665*47c6ae99SBarry Smith     Output Parameters:
666*47c6ae99SBarry Smith .   J  - matrix with the correct nonzero structure
667*47c6ae99SBarry Smith         (obviously without the correct Jacobian values)
668*47c6ae99SBarry Smith 
669*47c6ae99SBarry Smith     Level: advanced
670*47c6ae99SBarry Smith 
671*47c6ae99SBarry Smith     Notes: This properly preallocates the number of nonzeros in the sparse matrix so you
672*47c6ae99SBarry Smith        do not need to do it yourself.
673*47c6ae99SBarry Smith 
674*47c6ae99SBarry Smith        By default it also sets the nonzero structure and puts in the zero entries. To prevent setting
675*47c6ae99SBarry Smith        the nonzero pattern call DASetMatPreallocateOnly()
676*47c6ae99SBarry Smith 
677*47c6ae99SBarry Smith        When you call MatView() on this matrix it is displayed using the global natural ordering, NOT in the ordering used
678*47c6ae99SBarry Smith        internally by PETSc.
679*47c6ae99SBarry Smith 
680*47c6ae99SBarry Smith        In general it is easiest to use MatSetValuesStencil() or MatSetValuesLocal() to put values into the matrix because MatSetValues() requires
681*47c6ae99SBarry Smith        the indices for the global numbering for DAs which is complicated.
682*47c6ae99SBarry Smith 
683*47c6ae99SBarry Smith .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), DASetBlockFills(), DASetMatPreallocateOnly()
684*47c6ae99SBarry Smith 
685*47c6ae99SBarry Smith @*/
686*47c6ae99SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DAGetMatrix(DA da, const MatType mtype,Mat *J)
687*47c6ae99SBarry Smith {
688*47c6ae99SBarry Smith   PetscErrorCode ierr;
689*47c6ae99SBarry Smith   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
690*47c6ae99SBarry Smith   Mat            A;
691*47c6ae99SBarry Smith   MPI_Comm       comm;
692*47c6ae99SBarry Smith   const MatType  Atype;
693*47c6ae99SBarry Smith   void           (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL;
694*47c6ae99SBarry Smith   MatType        ttype[256];
695*47c6ae99SBarry Smith   PetscBool      flg;
696*47c6ae99SBarry Smith   PetscMPIInt    size;
697*47c6ae99SBarry Smith   DM_DA          *dd = (DM_DA*)da->data;
698*47c6ae99SBarry Smith 
699*47c6ae99SBarry Smith   PetscFunctionBegin;
700*47c6ae99SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES
701*47c6ae99SBarry Smith   ierr = MatInitializePackage(PETSC_NULL);CHKERRQ(ierr);
702*47c6ae99SBarry Smith #endif
703*47c6ae99SBarry Smith   ierr = PetscStrcpy((char*)ttype,mtype);CHKERRQ(ierr);
704*47c6ae99SBarry Smith   ierr = PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DA options","Mat");CHKERRQ(ierr);
705*47c6ae99SBarry Smith   ierr = PetscOptionsList("-da_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);CHKERRQ(ierr);
706*47c6ae99SBarry Smith   ierr = PetscOptionsEnd();
707*47c6ae99SBarry Smith 
708*47c6ae99SBarry Smith   /*
709*47c6ae99SBarry Smith                                   m
710*47c6ae99SBarry Smith           ------------------------------------------------------
711*47c6ae99SBarry Smith          |                                                     |
712*47c6ae99SBarry Smith          |                                                     |
713*47c6ae99SBarry Smith          |               ----------------------                |
714*47c6ae99SBarry Smith          |               |                    |                |
715*47c6ae99SBarry Smith       n  |           ny  |                    |                |
716*47c6ae99SBarry Smith          |               |                    |                |
717*47c6ae99SBarry Smith          |               .---------------------                |
718*47c6ae99SBarry Smith          |             (xs,ys)     nx                          |
719*47c6ae99SBarry Smith          |            .                                        |
720*47c6ae99SBarry Smith          |         (gxs,gys)                                   |
721*47c6ae99SBarry Smith          |                                                     |
722*47c6ae99SBarry Smith           -----------------------------------------------------
723*47c6ae99SBarry Smith   */
724*47c6ae99SBarry Smith 
725*47c6ae99SBarry Smith   /*
726*47c6ae99SBarry Smith          nc - number of components per grid point
727*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
728*47c6ae99SBarry Smith 
729*47c6ae99SBarry Smith   */
730*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0);CHKERRQ(ierr);
731*47c6ae99SBarry Smith   ierr = DAGetCorners(da,0,0,0,&nx,&ny,&nz);CHKERRQ(ierr);
732*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
733*47c6ae99SBarry Smith   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
734*47c6ae99SBarry Smith   ierr = MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);CHKERRQ(ierr);
735*47c6ae99SBarry Smith   ierr = MatSetType(A,(const MatType)ttype);CHKERRQ(ierr);
736*47c6ae99SBarry Smith   ierr = MatSetDA(A,da);CHKERRQ(ierr);
737*47c6ae99SBarry Smith   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
738*47c6ae99SBarry Smith   ierr = MatGetType(A,&Atype);CHKERRQ(ierr);
739*47c6ae99SBarry Smith   /*
740*47c6ae99SBarry Smith      We do not provide a getmatrix function in the DA operations because
741*47c6ae99SBarry Smith    the basic DA does not know about matrices. We think of DA as being more
742*47c6ae99SBarry Smith    more low-level than matrices. This is kind of cheating but, cause sometimes
743*47c6ae99SBarry Smith    we think of DA has higher level than matrices.
744*47c6ae99SBarry Smith 
745*47c6ae99SBarry Smith      We could switch based on Atype (or mtype), but we do not since the
746*47c6ae99SBarry Smith    specialized setting routines depend only the particular preallocation
747*47c6ae99SBarry Smith    details of the matrix, not the type itself.
748*47c6ae99SBarry Smith   */
749*47c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
750*47c6ae99SBarry Smith   if (!aij) {
751*47c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
752*47c6ae99SBarry Smith   }
753*47c6ae99SBarry Smith   if (!aij) {
754*47c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
755*47c6ae99SBarry Smith     if (!baij) {
756*47c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);CHKERRQ(ierr);
757*47c6ae99SBarry Smith     }
758*47c6ae99SBarry Smith     if (!baij){
759*47c6ae99SBarry Smith       ierr = PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
760*47c6ae99SBarry Smith       if (!sbaij) {
761*47c6ae99SBarry Smith         ierr = PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);CHKERRQ(ierr);
762*47c6ae99SBarry Smith       }
763*47c6ae99SBarry Smith     }
764*47c6ae99SBarry Smith   }
765*47c6ae99SBarry Smith   if (aij) {
766*47c6ae99SBarry Smith     if (dim == 1) {
767*47c6ae99SBarry Smith       ierr = DAGetMatrix1d_MPIAIJ(da,A);CHKERRQ(ierr);
768*47c6ae99SBarry Smith     } else if (dim == 2) {
769*47c6ae99SBarry Smith       if (dd->ofill) {
770*47c6ae99SBarry Smith         ierr = DAGetMatrix2d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
771*47c6ae99SBarry Smith       } else {
772*47c6ae99SBarry Smith         ierr = DAGetMatrix2d_MPIAIJ(da,A);CHKERRQ(ierr);
773*47c6ae99SBarry Smith       }
774*47c6ae99SBarry Smith     } else if (dim == 3) {
775*47c6ae99SBarry Smith       if (dd->ofill) {
776*47c6ae99SBarry Smith         ierr = DAGetMatrix3d_MPIAIJ_Fill(da,A);CHKERRQ(ierr);
777*47c6ae99SBarry Smith       } else {
778*47c6ae99SBarry Smith         ierr = DAGetMatrix3d_MPIAIJ(da,A);CHKERRQ(ierr);
779*47c6ae99SBarry Smith       }
780*47c6ae99SBarry Smith     }
781*47c6ae99SBarry Smith   } else if (baij) {
782*47c6ae99SBarry Smith     if (dim == 2) {
783*47c6ae99SBarry Smith       ierr = DAGetMatrix2d_MPIBAIJ(da,A);CHKERRQ(ierr);
784*47c6ae99SBarry Smith     } else if (dim == 3) {
785*47c6ae99SBarry Smith       ierr = DAGetMatrix3d_MPIBAIJ(da,A);CHKERRQ(ierr);
786*47c6ae99SBarry Smith     } else {
787*47c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \
788*47c6ae99SBarry Smith 	       "Send mail to petsc-maint@mcs.anl.gov for code",Atype,dim);
789*47c6ae99SBarry Smith     }
790*47c6ae99SBarry Smith   } else if (sbaij) {
791*47c6ae99SBarry Smith     if (dim == 2) {
792*47c6ae99SBarry Smith       ierr = DAGetMatrix2d_MPISBAIJ(da,A);CHKERRQ(ierr);
793*47c6ae99SBarry Smith     } else if (dim == 3) {
794*47c6ae99SBarry Smith       ierr = DAGetMatrix3d_MPISBAIJ(da,A);CHKERRQ(ierr);
795*47c6ae99SBarry Smith     } else {
796*47c6ae99SBarry Smith       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \
797*47c6ae99SBarry Smith 	       "Send mail to petsc-maint@mcs.anl.gov for code",Atype,dim);
798*47c6ae99SBarry Smith     }
799*47c6ae99SBarry Smith   }
800*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);CHKERRQ(ierr);
801*47c6ae99SBarry Smith   ierr = MatSetStencil(A,dim,dims,starts,dof);CHKERRQ(ierr);
802*47c6ae99SBarry Smith   ierr = PetscObjectCompose((PetscObject)A,"DA",(PetscObject)da);CHKERRQ(ierr);
803*47c6ae99SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
804*47c6ae99SBarry Smith   if (size > 1) {
805*47c6ae99SBarry Smith     /* change viewer to display matrix in natural ordering */
806*47c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);CHKERRQ(ierr);
807*47c6ae99SBarry Smith     ierr = MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);CHKERRQ(ierr);
808*47c6ae99SBarry Smith   }
809*47c6ae99SBarry Smith   *J = A;
810*47c6ae99SBarry Smith   PetscFunctionReturn(0);
811*47c6ae99SBarry Smith }
812*47c6ae99SBarry Smith 
813*47c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
814*47c6ae99SBarry Smith #undef __FUNCT__
815*47c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix2d_MPIAIJ"
816*47c6ae99SBarry Smith PetscErrorCode DAGetMatrix2d_MPIAIJ(DA da,Mat J)
817*47c6ae99SBarry Smith {
818*47c6ae99SBarry Smith   PetscErrorCode         ierr;
819*47c6ae99SBarry 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;
820*47c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
821*47c6ae99SBarry Smith   MPI_Comm               comm;
822*47c6ae99SBarry Smith   PetscScalar            *values;
823*47c6ae99SBarry Smith   DAPeriodicType         wrap;
824*47c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
825*47c6ae99SBarry Smith   DAStencilType          st;
826*47c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
827*47c6ae99SBarry Smith 
828*47c6ae99SBarry Smith   PetscFunctionBegin;
829*47c6ae99SBarry Smith   /*
830*47c6ae99SBarry Smith          nc - number of components per grid point
831*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
832*47c6ae99SBarry Smith 
833*47c6ae99SBarry Smith   */
834*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
835*47c6ae99SBarry Smith   col = 2*s + 1;
836*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
837*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
838*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
839*47c6ae99SBarry Smith 
840*47c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
841*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
842*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
843*47c6ae99SBarry Smith 
844*47c6ae99SBarry Smith   /* determine the matrix preallocation information */
845*47c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
846*47c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
847*47c6ae99SBarry Smith 
848*47c6ae99SBarry Smith     pstart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
849*47c6ae99SBarry Smith     pend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
850*47c6ae99SBarry Smith 
851*47c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
852*47c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
853*47c6ae99SBarry Smith 
854*47c6ae99SBarry Smith       lstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
855*47c6ae99SBarry Smith       lend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
856*47c6ae99SBarry Smith 
857*47c6ae99SBarry Smith       cnt  = 0;
858*47c6ae99SBarry Smith       for (k=0; k<nc; k++) {
859*47c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
860*47c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
861*47c6ae99SBarry Smith 	    if ((st == DA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
862*47c6ae99SBarry Smith 	      cols[cnt++]  = k + nc*(slot + gnx*l + p);
863*47c6ae99SBarry Smith 	    }
864*47c6ae99SBarry Smith 	  }
865*47c6ae99SBarry Smith 	}
866*47c6ae99SBarry Smith 	rows[k] = k + nc*(slot);
867*47c6ae99SBarry Smith       }
868*47c6ae99SBarry Smith       ierr = MatPreallocateSetLocal(ltog,nc,rows,cnt,cols,dnz,onz);CHKERRQ(ierr);
869*47c6ae99SBarry Smith     }
870*47c6ae99SBarry Smith   }
871*47c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
872*47c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
873*47c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
874*47c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
875*47c6ae99SBarry Smith 
876*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
877*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
878*47c6ae99SBarry Smith 
879*47c6ae99SBarry Smith   /*
880*47c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
881*47c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
882*47c6ae99SBarry Smith     PETSc ordering.
883*47c6ae99SBarry Smith   */
884*47c6ae99SBarry Smith   if (!dd->prealloc_only) {
885*47c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
886*47c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
887*47c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
888*47c6ae99SBarry Smith 
889*47c6ae99SBarry Smith       pstart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
890*47c6ae99SBarry Smith       pend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
891*47c6ae99SBarry Smith 
892*47c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
893*47c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
894*47c6ae99SBarry Smith 
895*47c6ae99SBarry Smith 	lstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
896*47c6ae99SBarry Smith 	lend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
897*47c6ae99SBarry Smith 
898*47c6ae99SBarry Smith 	cnt  = 0;
899*47c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
900*47c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
901*47c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
902*47c6ae99SBarry Smith 	      if ((st == DA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
903*47c6ae99SBarry Smith 		cols[cnt++]  = k + nc*(slot + gnx*l + p);
904*47c6ae99SBarry Smith 	      }
905*47c6ae99SBarry Smith 	    }
906*47c6ae99SBarry Smith 	  }
907*47c6ae99SBarry Smith 	  rows[k]      = k + nc*(slot);
908*47c6ae99SBarry Smith 	}
909*47c6ae99SBarry Smith 	ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
910*47c6ae99SBarry Smith       }
911*47c6ae99SBarry Smith     }
912*47c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
913*47c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
914*47c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
915*47c6ae99SBarry Smith   }
916*47c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
917*47c6ae99SBarry Smith   PetscFunctionReturn(0);
918*47c6ae99SBarry Smith }
919*47c6ae99SBarry Smith 
920*47c6ae99SBarry Smith #undef __FUNCT__
921*47c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix2d_MPIAIJ_Fill"
922*47c6ae99SBarry Smith PetscErrorCode DAGetMatrix2d_MPIAIJ_Fill(DA da,Mat J)
923*47c6ae99SBarry Smith {
924*47c6ae99SBarry Smith   PetscErrorCode         ierr;
925*47c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
926*47c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p;
927*47c6ae99SBarry Smith   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
928*47c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
929*47c6ae99SBarry Smith   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
930*47c6ae99SBarry Smith   MPI_Comm               comm;
931*47c6ae99SBarry Smith   PetscScalar            *values;
932*47c6ae99SBarry Smith   DAPeriodicType         wrap;
933*47c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
934*47c6ae99SBarry Smith   DAStencilType          st;
935*47c6ae99SBarry Smith 
936*47c6ae99SBarry Smith   PetscFunctionBegin;
937*47c6ae99SBarry Smith   /*
938*47c6ae99SBarry Smith          nc - number of components per grid point
939*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
940*47c6ae99SBarry Smith 
941*47c6ae99SBarry Smith   */
942*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
943*47c6ae99SBarry Smith   col = 2*s + 1;
944*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
945*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
946*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
947*47c6ae99SBarry Smith 
948*47c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
949*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
950*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
951*47c6ae99SBarry Smith 
952*47c6ae99SBarry Smith   /* determine the matrix preallocation information */
953*47c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);CHKERRQ(ierr);
954*47c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
955*47c6ae99SBarry Smith 
956*47c6ae99SBarry Smith     pstart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
957*47c6ae99SBarry Smith     pend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
958*47c6ae99SBarry Smith 
959*47c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
960*47c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
961*47c6ae99SBarry Smith 
962*47c6ae99SBarry Smith       lstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
963*47c6ae99SBarry Smith       lend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
964*47c6ae99SBarry Smith 
965*47c6ae99SBarry Smith       for (k=0; k<nc; k++) {
966*47c6ae99SBarry Smith         cnt  = 0;
967*47c6ae99SBarry Smith 	for (l=lstart; l<lend+1; l++) {
968*47c6ae99SBarry Smith 	  for (p=pstart; p<pend+1; p++) {
969*47c6ae99SBarry Smith             if (l || p) {
970*47c6ae99SBarry Smith 	      if ((st == DA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
971*47c6ae99SBarry Smith                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
972*47c6ae99SBarry Smith 		  cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
973*47c6ae99SBarry Smith 	      }
974*47c6ae99SBarry Smith             } else {
975*47c6ae99SBarry Smith 	      if (dfill) {
976*47c6ae99SBarry Smith 		for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
977*47c6ae99SBarry Smith 		  cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
978*47c6ae99SBarry Smith 	      } else {
979*47c6ae99SBarry Smith 		for (ifill_col=0; ifill_col<nc; ifill_col++)
980*47c6ae99SBarry Smith 		  cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
981*47c6ae99SBarry Smith 	      }
982*47c6ae99SBarry Smith             }
983*47c6ae99SBarry Smith 	  }
984*47c6ae99SBarry Smith 	}
985*47c6ae99SBarry Smith 	row = k + nc*(slot);
986*47c6ae99SBarry Smith         ierr = MatPreallocateSetLocal(ltog,1,&row,cnt,cols,dnz,onz);CHKERRQ(ierr);
987*47c6ae99SBarry Smith       }
988*47c6ae99SBarry Smith     }
989*47c6ae99SBarry Smith   }
990*47c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
991*47c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
992*47c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
993*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
994*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
995*47c6ae99SBarry Smith 
996*47c6ae99SBarry Smith   /*
997*47c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
998*47c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
999*47c6ae99SBarry Smith     PETSc ordering.
1000*47c6ae99SBarry Smith   */
1001*47c6ae99SBarry Smith   if (!dd->prealloc_only) {
1002*47c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1003*47c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1004*47c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1005*47c6ae99SBarry Smith 
1006*47c6ae99SBarry Smith       pstart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1007*47c6ae99SBarry Smith       pend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1008*47c6ae99SBarry Smith 
1009*47c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1010*47c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
1011*47c6ae99SBarry Smith 
1012*47c6ae99SBarry Smith 	lstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1013*47c6ae99SBarry Smith 	lend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1014*47c6ae99SBarry Smith 
1015*47c6ae99SBarry Smith 	for (k=0; k<nc; k++) {
1016*47c6ae99SBarry Smith 	  cnt  = 0;
1017*47c6ae99SBarry Smith 	  for (l=lstart; l<lend+1; l++) {
1018*47c6ae99SBarry Smith 	    for (p=pstart; p<pend+1; p++) {
1019*47c6ae99SBarry Smith 	      if (l || p) {
1020*47c6ae99SBarry Smith 		if ((st == DA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
1021*47c6ae99SBarry Smith 		  for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
1022*47c6ae99SBarry Smith 		    cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
1023*47c6ae99SBarry Smith 		}
1024*47c6ae99SBarry Smith 	      } else {
1025*47c6ae99SBarry Smith 		if (dfill) {
1026*47c6ae99SBarry Smith 		  for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
1027*47c6ae99SBarry Smith 		    cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
1028*47c6ae99SBarry Smith 		} else {
1029*47c6ae99SBarry Smith 		  for (ifill_col=0; ifill_col<nc; ifill_col++)
1030*47c6ae99SBarry Smith 		    cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
1031*47c6ae99SBarry Smith 		}
1032*47c6ae99SBarry Smith 	      }
1033*47c6ae99SBarry Smith 	    }
1034*47c6ae99SBarry Smith 	  }
1035*47c6ae99SBarry Smith 	  row  = k + nc*(slot);
1036*47c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1037*47c6ae99SBarry Smith 	}
1038*47c6ae99SBarry Smith       }
1039*47c6ae99SBarry Smith     }
1040*47c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1041*47c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1042*47c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1043*47c6ae99SBarry Smith   }
1044*47c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
1045*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1046*47c6ae99SBarry Smith }
1047*47c6ae99SBarry Smith 
1048*47c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
1049*47c6ae99SBarry Smith 
1050*47c6ae99SBarry Smith #undef __FUNCT__
1051*47c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix3d_MPIAIJ"
1052*47c6ae99SBarry Smith PetscErrorCode DAGetMatrix3d_MPIAIJ(DA da,Mat J)
1053*47c6ae99SBarry Smith {
1054*47c6ae99SBarry Smith   PetscErrorCode         ierr;
1055*47c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1056*47c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL;
1057*47c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
1058*47c6ae99SBarry Smith   MPI_Comm               comm;
1059*47c6ae99SBarry Smith   PetscScalar            *values;
1060*47c6ae99SBarry Smith   DAPeriodicType         wrap;
1061*47c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1062*47c6ae99SBarry Smith   DAStencilType          st;
1063*47c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1064*47c6ae99SBarry Smith 
1065*47c6ae99SBarry Smith   PetscFunctionBegin;
1066*47c6ae99SBarry Smith   /*
1067*47c6ae99SBarry Smith          nc - number of components per grid point
1068*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
1069*47c6ae99SBarry Smith 
1070*47c6ae99SBarry Smith   */
1071*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
1072*47c6ae99SBarry Smith   col    = 2*s + 1;
1073*47c6ae99SBarry Smith 
1074*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1075*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
1076*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
1077*47c6ae99SBarry Smith 
1078*47c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
1079*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1080*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
1081*47c6ae99SBarry Smith 
1082*47c6ae99SBarry Smith   /* determine the matrix preallocation information */
1083*47c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
1084*47c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1085*47c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1086*47c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1087*47c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1088*47c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1089*47c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1090*47c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1091*47c6ae99SBarry Smith 	kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
1092*47c6ae99SBarry Smith 	kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
1093*47c6ae99SBarry Smith 
1094*47c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1095*47c6ae99SBarry Smith 
1096*47c6ae99SBarry Smith 	cnt  = 0;
1097*47c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
1098*47c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
1099*47c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
1100*47c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
1101*47c6ae99SBarry Smith 		if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1102*47c6ae99SBarry Smith 		  cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1103*47c6ae99SBarry Smith 		}
1104*47c6ae99SBarry Smith 	      }
1105*47c6ae99SBarry Smith 	    }
1106*47c6ae99SBarry Smith 	  }
1107*47c6ae99SBarry Smith 	  rows[l] = l + nc*(slot);
1108*47c6ae99SBarry Smith 	}
1109*47c6ae99SBarry Smith 	ierr = MatPreallocateSetLocal(ltog,nc,rows,cnt,cols,dnz,onz);CHKERRQ(ierr);
1110*47c6ae99SBarry Smith       }
1111*47c6ae99SBarry Smith     }
1112*47c6ae99SBarry Smith   }
1113*47c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
1114*47c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
1115*47c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1116*47c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
1117*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
1118*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
1119*47c6ae99SBarry Smith 
1120*47c6ae99SBarry Smith   /*
1121*47c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1122*47c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1123*47c6ae99SBarry Smith     PETSc ordering.
1124*47c6ae99SBarry Smith   */
1125*47c6ae99SBarry Smith   if (!dd->prealloc_only) {
1126*47c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1127*47c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1128*47c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1129*47c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1130*47c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1131*47c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1132*47c6ae99SBarry Smith 	jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1133*47c6ae99SBarry Smith 	jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1134*47c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
1135*47c6ae99SBarry Smith 	  kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
1136*47c6ae99SBarry Smith 	  kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
1137*47c6ae99SBarry Smith 
1138*47c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1139*47c6ae99SBarry Smith 
1140*47c6ae99SBarry Smith 	  cnt  = 0;
1141*47c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
1142*47c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
1143*47c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
1144*47c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
1145*47c6ae99SBarry Smith 		  if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1146*47c6ae99SBarry Smith 		    cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1147*47c6ae99SBarry Smith 		  }
1148*47c6ae99SBarry Smith 		}
1149*47c6ae99SBarry Smith 	      }
1150*47c6ae99SBarry Smith 	    }
1151*47c6ae99SBarry Smith 	    rows[l]      = l + nc*(slot);
1152*47c6ae99SBarry Smith 	  }
1153*47c6ae99SBarry Smith 	  ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1154*47c6ae99SBarry Smith 	}
1155*47c6ae99SBarry Smith       }
1156*47c6ae99SBarry Smith     }
1157*47c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1158*47c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1159*47c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1160*47c6ae99SBarry Smith   }
1161*47c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
1162*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1163*47c6ae99SBarry Smith }
1164*47c6ae99SBarry Smith 
1165*47c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
1166*47c6ae99SBarry Smith 
1167*47c6ae99SBarry Smith #undef __FUNCT__
1168*47c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix1d_MPIAIJ"
1169*47c6ae99SBarry Smith PetscErrorCode DAGetMatrix1d_MPIAIJ(DA da,Mat J)
1170*47c6ae99SBarry Smith {
1171*47c6ae99SBarry Smith   PetscErrorCode         ierr;
1172*47c6ae99SBarry Smith   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
1173*47c6ae99SBarry Smith   PetscInt               m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l;
1174*47c6ae99SBarry Smith   PetscInt               istart,iend;
1175*47c6ae99SBarry Smith   PetscScalar            *values;
1176*47c6ae99SBarry Smith   DAPeriodicType         wrap;
1177*47c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1178*47c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1179*47c6ae99SBarry Smith 
1180*47c6ae99SBarry Smith   PetscFunctionBegin;
1181*47c6ae99SBarry Smith   /*
1182*47c6ae99SBarry Smith          nc - number of components per grid point
1183*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
1184*47c6ae99SBarry Smith 
1185*47c6ae99SBarry Smith   */
1186*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&wrap,0);CHKERRQ(ierr);
1187*47c6ae99SBarry Smith   col    = 2*s + 1;
1188*47c6ae99SBarry Smith 
1189*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,0,0,&nx,0,0);CHKERRQ(ierr);
1190*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);CHKERRQ(ierr);
1191*47c6ae99SBarry Smith 
1192*47c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,col*nc,0);CHKERRQ(ierr);
1193*47c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);CHKERRQ(ierr);
1194*47c6ae99SBarry Smith   ierr = MatSetBlockSize(J,nc);CHKERRQ(ierr);
1195*47c6ae99SBarry Smith   ierr = PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);CHKERRQ(ierr);
1196*47c6ae99SBarry Smith 
1197*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1198*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
1199*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
1200*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
1201*47c6ae99SBarry Smith 
1202*47c6ae99SBarry Smith   /*
1203*47c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1204*47c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1205*47c6ae99SBarry Smith     PETSc ordering.
1206*47c6ae99SBarry Smith   */
1207*47c6ae99SBarry Smith   if (!dd->prealloc_only) {
1208*47c6ae99SBarry Smith     ierr = PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1209*47c6ae99SBarry Smith     ierr = PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1210*47c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1211*47c6ae99SBarry Smith       istart = PetscMax(-s,gxs - i);
1212*47c6ae99SBarry Smith       iend   = PetscMin(s,gxs + gnx - i - 1);
1213*47c6ae99SBarry Smith       slot   = i - gxs;
1214*47c6ae99SBarry Smith 
1215*47c6ae99SBarry Smith       cnt  = 0;
1216*47c6ae99SBarry Smith       for (l=0; l<nc; l++) {
1217*47c6ae99SBarry Smith 	for (i1=istart; i1<iend+1; i1++) {
1218*47c6ae99SBarry Smith 	  cols[cnt++] = l + nc*(slot + i1);
1219*47c6ae99SBarry Smith 	}
1220*47c6ae99SBarry Smith 	rows[l]      = l + nc*(slot);
1221*47c6ae99SBarry Smith       }
1222*47c6ae99SBarry Smith       ierr = MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1223*47c6ae99SBarry Smith     }
1224*47c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1225*47c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1226*47c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1227*47c6ae99SBarry Smith   }
1228*47c6ae99SBarry Smith   ierr = PetscFree2(rows,cols);CHKERRQ(ierr);
1229*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1230*47c6ae99SBarry Smith }
1231*47c6ae99SBarry Smith 
1232*47c6ae99SBarry Smith #undef __FUNCT__
1233*47c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix2d_MPIBAIJ"
1234*47c6ae99SBarry Smith PetscErrorCode DAGetMatrix2d_MPIBAIJ(DA da,Mat J)
1235*47c6ae99SBarry Smith {
1236*47c6ae99SBarry Smith   PetscErrorCode         ierr;
1237*47c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1238*47c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
1239*47c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
1240*47c6ae99SBarry Smith   MPI_Comm               comm;
1241*47c6ae99SBarry Smith   PetscScalar            *values;
1242*47c6ae99SBarry Smith   DAPeriodicType         wrap;
1243*47c6ae99SBarry Smith   DAStencilType          st;
1244*47c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1245*47c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1246*47c6ae99SBarry Smith 
1247*47c6ae99SBarry Smith   PetscFunctionBegin;
1248*47c6ae99SBarry Smith   /*
1249*47c6ae99SBarry Smith      nc - number of components per grid point
1250*47c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
1251*47c6ae99SBarry Smith   */
1252*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
1253*47c6ae99SBarry Smith   col = 2*s + 1;
1254*47c6ae99SBarry Smith 
1255*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1256*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
1257*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
1258*47c6ae99SBarry Smith 
1259*47c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1260*47c6ae99SBarry Smith 
1261*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1262*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
1263*47c6ae99SBarry Smith 
1264*47c6ae99SBarry Smith   /* determine the matrix preallocation information */
1265*47c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
1266*47c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1267*47c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1268*47c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1269*47c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1270*47c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1271*47c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1272*47c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
1273*47c6ae99SBarry Smith 
1274*47c6ae99SBarry Smith       /* Find block columns in block row */
1275*47c6ae99SBarry Smith       cnt  = 0;
1276*47c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
1277*47c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1278*47c6ae99SBarry Smith           if (st == DA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
1279*47c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
1280*47c6ae99SBarry Smith           }
1281*47c6ae99SBarry Smith         }
1282*47c6ae99SBarry Smith       }
1283*47c6ae99SBarry Smith       ierr = MatPreallocateSetLocal(ltogb,1,&slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
1284*47c6ae99SBarry Smith     }
1285*47c6ae99SBarry Smith   }
1286*47c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
1287*47c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
1288*47c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1289*47c6ae99SBarry Smith 
1290*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
1291*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
1292*47c6ae99SBarry Smith 
1293*47c6ae99SBarry Smith   /*
1294*47c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1295*47c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1296*47c6ae99SBarry Smith     PETSc ordering.
1297*47c6ae99SBarry Smith   */
1298*47c6ae99SBarry Smith   if (!dd->prealloc_only) {
1299*47c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1300*47c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1301*47c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1302*47c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1303*47c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1304*47c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1305*47c6ae99SBarry Smith 	jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1306*47c6ae99SBarry Smith 	jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1307*47c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys);
1308*47c6ae99SBarry Smith 	cnt  = 0;
1309*47c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
1310*47c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1311*47c6ae99SBarry Smith             if (st == DA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
1312*47c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
1313*47c6ae99SBarry Smith             }
1314*47c6ae99SBarry Smith           }
1315*47c6ae99SBarry Smith         }
1316*47c6ae99SBarry Smith 	ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1317*47c6ae99SBarry Smith       }
1318*47c6ae99SBarry Smith     }
1319*47c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1320*47c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1321*47c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1322*47c6ae99SBarry Smith   }
1323*47c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
1324*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1325*47c6ae99SBarry Smith }
1326*47c6ae99SBarry Smith 
1327*47c6ae99SBarry Smith #undef __FUNCT__
1328*47c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix3d_MPIBAIJ"
1329*47c6ae99SBarry Smith PetscErrorCode DAGetMatrix3d_MPIBAIJ(DA da,Mat J)
1330*47c6ae99SBarry Smith {
1331*47c6ae99SBarry Smith   PetscErrorCode         ierr;
1332*47c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1333*47c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
1334*47c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
1335*47c6ae99SBarry Smith   MPI_Comm               comm;
1336*47c6ae99SBarry Smith   PetscScalar            *values;
1337*47c6ae99SBarry Smith   DAPeriodicType         wrap;
1338*47c6ae99SBarry Smith   DAStencilType          st;
1339*47c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1340*47c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1341*47c6ae99SBarry Smith 
1342*47c6ae99SBarry Smith   PetscFunctionBegin;
1343*47c6ae99SBarry Smith   /*
1344*47c6ae99SBarry Smith          nc - number of components per grid point
1345*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
1346*47c6ae99SBarry Smith 
1347*47c6ae99SBarry Smith   */
1348*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
1349*47c6ae99SBarry Smith   col    = 2*s + 1;
1350*47c6ae99SBarry Smith 
1351*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1352*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
1353*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
1354*47c6ae99SBarry Smith 
1355*47c6ae99SBarry Smith   ierr  = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1356*47c6ae99SBarry Smith 
1357*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1358*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
1359*47c6ae99SBarry Smith 
1360*47c6ae99SBarry Smith   /* determine the matrix preallocation information */
1361*47c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
1362*47c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1363*47c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1364*47c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1365*47c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1366*47c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1367*47c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1368*47c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1369*47c6ae99SBarry Smith 	kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
1370*47c6ae99SBarry Smith 	kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
1371*47c6ae99SBarry Smith 
1372*47c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1373*47c6ae99SBarry Smith 
1374*47c6ae99SBarry Smith 	/* Find block columns in block row */
1375*47c6ae99SBarry Smith 	cnt  = 0;
1376*47c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
1377*47c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1378*47c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1379*47c6ae99SBarry Smith               if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1380*47c6ae99SBarry Smith 		cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
1381*47c6ae99SBarry Smith 	      }
1382*47c6ae99SBarry Smith 	    }
1383*47c6ae99SBarry Smith 	  }
1384*47c6ae99SBarry Smith 	}
1385*47c6ae99SBarry Smith 	ierr = MatPreallocateSetLocal(ltogb,1,&slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
1386*47c6ae99SBarry Smith       }
1387*47c6ae99SBarry Smith     }
1388*47c6ae99SBarry Smith   }
1389*47c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
1390*47c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
1391*47c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1392*47c6ae99SBarry Smith 
1393*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
1394*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
1395*47c6ae99SBarry Smith 
1396*47c6ae99SBarry Smith   /*
1397*47c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1398*47c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1399*47c6ae99SBarry Smith     PETSc ordering.
1400*47c6ae99SBarry Smith   */
1401*47c6ae99SBarry Smith   if (!dd->prealloc_only) {
1402*47c6ae99SBarry Smith     ierr  = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1403*47c6ae99SBarry Smith     ierr  = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1404*47c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1405*47c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1406*47c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1407*47c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1408*47c6ae99SBarry Smith 	jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1409*47c6ae99SBarry Smith 	jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1410*47c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
1411*47c6ae99SBarry Smith 	  kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
1412*47c6ae99SBarry Smith 	  kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
1413*47c6ae99SBarry Smith 
1414*47c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1415*47c6ae99SBarry Smith 
1416*47c6ae99SBarry Smith 	  cnt  = 0;
1417*47c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
1418*47c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
1419*47c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1420*47c6ae99SBarry Smith                 if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1421*47c6ae99SBarry Smith                   cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
1422*47c6ae99SBarry Smith                 }
1423*47c6ae99SBarry Smith               }
1424*47c6ae99SBarry Smith             }
1425*47c6ae99SBarry Smith           }
1426*47c6ae99SBarry Smith 	  ierr = MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1427*47c6ae99SBarry Smith 	}
1428*47c6ae99SBarry Smith       }
1429*47c6ae99SBarry Smith     }
1430*47c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1431*47c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1432*47c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1433*47c6ae99SBarry Smith   }
1434*47c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
1435*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1436*47c6ae99SBarry Smith }
1437*47c6ae99SBarry Smith 
1438*47c6ae99SBarry Smith #undef __FUNCT__
1439*47c6ae99SBarry Smith #define __FUNCT__ "L2GFilterUpperTriangular"
1440*47c6ae99SBarry Smith /*
1441*47c6ae99SBarry Smith   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
1442*47c6ae99SBarry Smith   identify in the local ordering with periodic domain.
1443*47c6ae99SBarry Smith */
1444*47c6ae99SBarry Smith static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
1445*47c6ae99SBarry Smith {
1446*47c6ae99SBarry Smith   PetscErrorCode ierr;
1447*47c6ae99SBarry Smith   PetscInt i,n;
1448*47c6ae99SBarry Smith 
1449*47c6ae99SBarry Smith   PetscFunctionBegin;
1450*47c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,1,row,row);CHKERRQ(ierr);
1451*47c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingApply(ltog,*cnt,col,col);CHKERRQ(ierr);
1452*47c6ae99SBarry Smith   for (i=0,n=0; i<*cnt; i++) {
1453*47c6ae99SBarry Smith     if (col[i] >= *row) col[n++] = col[i];
1454*47c6ae99SBarry Smith   }
1455*47c6ae99SBarry Smith   *cnt = n;
1456*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1457*47c6ae99SBarry Smith }
1458*47c6ae99SBarry Smith 
1459*47c6ae99SBarry Smith #undef __FUNCT__
1460*47c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix2d_MPISBAIJ"
1461*47c6ae99SBarry Smith PetscErrorCode DAGetMatrix2d_MPISBAIJ(DA da,Mat J)
1462*47c6ae99SBarry Smith {
1463*47c6ae99SBarry Smith   PetscErrorCode         ierr;
1464*47c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1465*47c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
1466*47c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,ii,jj;
1467*47c6ae99SBarry Smith   MPI_Comm               comm;
1468*47c6ae99SBarry Smith   PetscScalar            *values;
1469*47c6ae99SBarry Smith   DAPeriodicType         wrap;
1470*47c6ae99SBarry Smith   DAStencilType          st;
1471*47c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1472*47c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1473*47c6ae99SBarry Smith 
1474*47c6ae99SBarry Smith   PetscFunctionBegin;
1475*47c6ae99SBarry Smith   /*
1476*47c6ae99SBarry Smith      nc - number of components per grid point
1477*47c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
1478*47c6ae99SBarry Smith   */
1479*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
1480*47c6ae99SBarry Smith   col = 2*s + 1;
1481*47c6ae99SBarry Smith 
1482*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);CHKERRQ(ierr);
1483*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);CHKERRQ(ierr);
1484*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
1485*47c6ae99SBarry Smith 
1486*47c6ae99SBarry Smith   ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1487*47c6ae99SBarry Smith 
1488*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1489*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
1490*47c6ae99SBarry Smith 
1491*47c6ae99SBarry Smith   /* determine the matrix preallocation information */
1492*47c6ae99SBarry Smith   ierr = MatPreallocateSymmetricInitialize(comm,nx*ny,nx*ny,dnz,onz);CHKERRQ(ierr);
1493*47c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1494*47c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1495*47c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1496*47c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1497*47c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1498*47c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1499*47c6ae99SBarry Smith       slot = i - gxs + gnx*(j - gys);
1500*47c6ae99SBarry Smith 
1501*47c6ae99SBarry Smith       /* Find block columns in block row */
1502*47c6ae99SBarry Smith       cnt  = 0;
1503*47c6ae99SBarry Smith       for (ii=istart; ii<iend+1; ii++) {
1504*47c6ae99SBarry Smith         for (jj=jstart; jj<jend+1; jj++) {
1505*47c6ae99SBarry Smith           if (st == DA_STENCIL_BOX || !ii || !jj) {
1506*47c6ae99SBarry Smith             cols[cnt++]  = slot + ii + gnx*jj;
1507*47c6ae99SBarry Smith           }
1508*47c6ae99SBarry Smith         }
1509*47c6ae99SBarry Smith       }
1510*47c6ae99SBarry Smith       ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
1511*47c6ae99SBarry Smith       ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
1512*47c6ae99SBarry Smith     }
1513*47c6ae99SBarry Smith   }
1514*47c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
1515*47c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
1516*47c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1517*47c6ae99SBarry Smith 
1518*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
1519*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
1520*47c6ae99SBarry Smith 
1521*47c6ae99SBarry Smith   /*
1522*47c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1523*47c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1524*47c6ae99SBarry Smith     PETSc ordering.
1525*47c6ae99SBarry Smith   */
1526*47c6ae99SBarry Smith   if (!dd->prealloc_only) {
1527*47c6ae99SBarry Smith     ierr = PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1528*47c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1529*47c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1530*47c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1531*47c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1532*47c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1533*47c6ae99SBarry Smith         jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1534*47c6ae99SBarry Smith         jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1535*47c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys);
1536*47c6ae99SBarry Smith 
1537*47c6ae99SBarry Smith         /* Find block columns in block row */
1538*47c6ae99SBarry Smith         cnt  = 0;
1539*47c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
1540*47c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1541*47c6ae99SBarry Smith             if (st == DA_STENCIL_BOX || !ii || !jj) {
1542*47c6ae99SBarry Smith               cols[cnt++]  = slot + ii + gnx*jj;
1543*47c6ae99SBarry Smith             }
1544*47c6ae99SBarry Smith           }
1545*47c6ae99SBarry Smith         }
1546*47c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
1547*47c6ae99SBarry Smith 	ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1548*47c6ae99SBarry Smith       }
1549*47c6ae99SBarry Smith     }
1550*47c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1551*47c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1552*47c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1553*47c6ae99SBarry Smith   }
1554*47c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
1555*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1556*47c6ae99SBarry Smith }
1557*47c6ae99SBarry Smith 
1558*47c6ae99SBarry Smith #undef __FUNCT__
1559*47c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix3d_MPISBAIJ"
1560*47c6ae99SBarry Smith PetscErrorCode DAGetMatrix3d_MPISBAIJ(DA da,Mat J)
1561*47c6ae99SBarry Smith {
1562*47c6ae99SBarry Smith   PetscErrorCode         ierr;
1563*47c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1564*47c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
1565*47c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
1566*47c6ae99SBarry Smith   MPI_Comm               comm;
1567*47c6ae99SBarry Smith   PetscScalar            *values;
1568*47c6ae99SBarry Smith   DAPeriodicType         wrap;
1569*47c6ae99SBarry Smith   DAStencilType          st;
1570*47c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1571*47c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1572*47c6ae99SBarry Smith 
1573*47c6ae99SBarry Smith   PetscFunctionBegin;
1574*47c6ae99SBarry Smith   /*
1575*47c6ae99SBarry Smith      nc - number of components per grid point
1576*47c6ae99SBarry Smith      col - number of colors needed in one direction for single component problem
1577*47c6ae99SBarry Smith   */
1578*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
1579*47c6ae99SBarry Smith   col = 2*s + 1;
1580*47c6ae99SBarry Smith 
1581*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1582*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
1583*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
1584*47c6ae99SBarry Smith 
1585*47c6ae99SBarry Smith   /* create the matrix */
1586*47c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1587*47c6ae99SBarry Smith 
1588*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1589*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
1590*47c6ae99SBarry Smith 
1591*47c6ae99SBarry Smith   /* determine the matrix preallocation information */
1592*47c6ae99SBarry Smith   ierr = MatPreallocateSymmetricInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);CHKERRQ(ierr);
1593*47c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1594*47c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1595*47c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1596*47c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1597*47c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1598*47c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1599*47c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1600*47c6ae99SBarry Smith         kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
1601*47c6ae99SBarry Smith 	kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
1602*47c6ae99SBarry Smith 
1603*47c6ae99SBarry Smith 	slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1604*47c6ae99SBarry Smith 
1605*47c6ae99SBarry Smith 	/* Find block columns in block row */
1606*47c6ae99SBarry Smith 	cnt  = 0;
1607*47c6ae99SBarry Smith         for (ii=istart; ii<iend+1; ii++) {
1608*47c6ae99SBarry Smith           for (jj=jstart; jj<jend+1; jj++) {
1609*47c6ae99SBarry Smith             for (kk=kstart; kk<kend+1; kk++) {
1610*47c6ae99SBarry Smith               if ((st == DA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
1611*47c6ae99SBarry Smith                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
1612*47c6ae99SBarry Smith               }
1613*47c6ae99SBarry Smith             }
1614*47c6ae99SBarry Smith           }
1615*47c6ae99SBarry Smith         }
1616*47c6ae99SBarry Smith         ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
1617*47c6ae99SBarry Smith         ierr = MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);CHKERRQ(ierr);
1618*47c6ae99SBarry Smith       }
1619*47c6ae99SBarry Smith     }
1620*47c6ae99SBarry Smith   }
1621*47c6ae99SBarry Smith   ierr = MatSeqSBAIJSetPreallocation(J,nc,0,dnz);CHKERRQ(ierr);
1622*47c6ae99SBarry Smith   ierr = MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);CHKERRQ(ierr);
1623*47c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1624*47c6ae99SBarry Smith 
1625*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
1626*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
1627*47c6ae99SBarry Smith 
1628*47c6ae99SBarry Smith   /*
1629*47c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1630*47c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1631*47c6ae99SBarry Smith     PETSc ordering.
1632*47c6ae99SBarry Smith   */
1633*47c6ae99SBarry Smith   if (!dd->prealloc_only) {
1634*47c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1635*47c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1636*47c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1637*47c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1638*47c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1639*47c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1640*47c6ae99SBarry Smith         jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1641*47c6ae99SBarry Smith 	jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1642*47c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
1643*47c6ae99SBarry Smith           kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
1644*47c6ae99SBarry Smith 	  kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
1645*47c6ae99SBarry Smith 
1646*47c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1647*47c6ae99SBarry Smith 
1648*47c6ae99SBarry Smith 	  cnt  = 0;
1649*47c6ae99SBarry Smith           for (ii=istart; ii<iend+1; ii++) {
1650*47c6ae99SBarry Smith             for (jj=jstart; jj<jend+1; jj++) {
1651*47c6ae99SBarry Smith               for (kk=kstart; kk<kend+1; kk++) {
1652*47c6ae99SBarry Smith                 if ((st == DA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
1653*47c6ae99SBarry Smith 		  cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
1654*47c6ae99SBarry Smith 		}
1655*47c6ae99SBarry Smith 	      }
1656*47c6ae99SBarry Smith 	    }
1657*47c6ae99SBarry Smith 	  }
1658*47c6ae99SBarry Smith           ierr = L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);CHKERRQ(ierr);
1659*47c6ae99SBarry Smith           ierr = MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1660*47c6ae99SBarry Smith 	}
1661*47c6ae99SBarry Smith       }
1662*47c6ae99SBarry Smith     }
1663*47c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1664*47c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1665*47c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1666*47c6ae99SBarry Smith   }
1667*47c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
1668*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1669*47c6ae99SBarry Smith }
1670*47c6ae99SBarry Smith 
1671*47c6ae99SBarry Smith /* ---------------------------------------------------------------------------------*/
1672*47c6ae99SBarry Smith 
1673*47c6ae99SBarry Smith #undef __FUNCT__
1674*47c6ae99SBarry Smith #define __FUNCT__ "DAGetMatrix3d_MPIAIJ_Fill"
1675*47c6ae99SBarry Smith PetscErrorCode DAGetMatrix3d_MPIAIJ_Fill(DA da,Mat J)
1676*47c6ae99SBarry Smith {
1677*47c6ae99SBarry Smith   PetscErrorCode         ierr;
1678*47c6ae99SBarry Smith   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1679*47c6ae99SBarry Smith   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz;
1680*47c6ae99SBarry Smith   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
1681*47c6ae99SBarry Smith   DM_DA                  *dd = (DM_DA*)da->data;
1682*47c6ae99SBarry Smith   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
1683*47c6ae99SBarry Smith   MPI_Comm               comm;
1684*47c6ae99SBarry Smith   PetscScalar            *values;
1685*47c6ae99SBarry Smith   DAPeriodicType         wrap;
1686*47c6ae99SBarry Smith   ISLocalToGlobalMapping ltog,ltogb;
1687*47c6ae99SBarry Smith   DAStencilType          st;
1688*47c6ae99SBarry Smith 
1689*47c6ae99SBarry Smith   PetscFunctionBegin;
1690*47c6ae99SBarry Smith   /*
1691*47c6ae99SBarry Smith          nc - number of components per grid point
1692*47c6ae99SBarry Smith          col - number of colors needed in one direction for single component problem
1693*47c6ae99SBarry Smith 
1694*47c6ae99SBarry Smith   */
1695*47c6ae99SBarry Smith   ierr = DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);CHKERRQ(ierr);
1696*47c6ae99SBarry Smith   col    = 2*s + 1;
1697*47c6ae99SBarry Smith   if (DAXPeriodic(wrap) && (m % col)){
1698*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
1699*47c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1700*47c6ae99SBarry Smith   }
1701*47c6ae99SBarry Smith   if (DAYPeriodic(wrap) && (n % col)){
1702*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
1703*47c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1704*47c6ae99SBarry Smith   }
1705*47c6ae99SBarry Smith   if (DAZPeriodic(wrap) && (p % col)){
1706*47c6ae99SBarry Smith     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
1707*47c6ae99SBarry Smith                  by 2*stencil_width + 1\n");
1708*47c6ae99SBarry Smith   }
1709*47c6ae99SBarry Smith 
1710*47c6ae99SBarry Smith   ierr = DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);CHKERRQ(ierr);
1711*47c6ae99SBarry Smith   ierr = DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);CHKERRQ(ierr);
1712*47c6ae99SBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
1713*47c6ae99SBarry Smith 
1714*47c6ae99SBarry Smith   ierr = PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);CHKERRQ(ierr);
1715*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMapping(da,&ltog);CHKERRQ(ierr);
1716*47c6ae99SBarry Smith   ierr = DAGetISLocalToGlobalMappingBlck(da,&ltogb);CHKERRQ(ierr);
1717*47c6ae99SBarry Smith 
1718*47c6ae99SBarry Smith   /* determine the matrix preallocation information */
1719*47c6ae99SBarry Smith   ierr = MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);CHKERRQ(ierr);
1720*47c6ae99SBarry Smith 
1721*47c6ae99SBarry Smith 
1722*47c6ae99SBarry Smith   for (i=xs; i<xs+nx; i++) {
1723*47c6ae99SBarry Smith     istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1724*47c6ae99SBarry Smith     iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1725*47c6ae99SBarry Smith     for (j=ys; j<ys+ny; j++) {
1726*47c6ae99SBarry Smith       jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1727*47c6ae99SBarry Smith       jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1728*47c6ae99SBarry Smith       for (k=zs; k<zs+nz; k++) {
1729*47c6ae99SBarry Smith         kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
1730*47c6ae99SBarry Smith         kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
1731*47c6ae99SBarry Smith 
1732*47c6ae99SBarry Smith         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1733*47c6ae99SBarry Smith 
1734*47c6ae99SBarry Smith 	for (l=0; l<nc; l++) {
1735*47c6ae99SBarry Smith 	  cnt  = 0;
1736*47c6ae99SBarry Smith 	  for (ii=istart; ii<iend+1; ii++) {
1737*47c6ae99SBarry Smith 	    for (jj=jstart; jj<jend+1; jj++) {
1738*47c6ae99SBarry Smith 	      for (kk=kstart; kk<kend+1; kk++) {
1739*47c6ae99SBarry Smith 		if (ii || jj || kk) {
1740*47c6ae99SBarry Smith 		  if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1741*47c6ae99SBarry Smith 		    for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
1742*47c6ae99SBarry Smith 		      cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1743*47c6ae99SBarry Smith 		  }
1744*47c6ae99SBarry Smith 		} else {
1745*47c6ae99SBarry Smith 		  if (dfill) {
1746*47c6ae99SBarry Smith 		    for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
1747*47c6ae99SBarry Smith 		      cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1748*47c6ae99SBarry Smith 		  } else {
1749*47c6ae99SBarry Smith 		    for (ifill_col=0; ifill_col<nc; ifill_col++)
1750*47c6ae99SBarry Smith 		      cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1751*47c6ae99SBarry Smith 		  }
1752*47c6ae99SBarry Smith 		}
1753*47c6ae99SBarry Smith 	      }
1754*47c6ae99SBarry Smith 	    }
1755*47c6ae99SBarry Smith 	  }
1756*47c6ae99SBarry Smith 	  row  = l + nc*(slot);
1757*47c6ae99SBarry Smith 	  ierr = MatPreallocateSetLocal(ltog,1,&row,cnt,cols,dnz,onz);CHKERRQ(ierr);
1758*47c6ae99SBarry Smith 	}
1759*47c6ae99SBarry Smith       }
1760*47c6ae99SBarry Smith     }
1761*47c6ae99SBarry Smith   }
1762*47c6ae99SBarry Smith   ierr = MatSeqAIJSetPreallocation(J,0,dnz);CHKERRQ(ierr);
1763*47c6ae99SBarry Smith   ierr = MatMPIAIJSetPreallocation(J,0,dnz,0,onz);CHKERRQ(ierr);
1764*47c6ae99SBarry Smith   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
1765*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMapping(J,ltog);CHKERRQ(ierr);
1766*47c6ae99SBarry Smith   ierr = MatSetLocalToGlobalMappingBlock(J,ltogb);CHKERRQ(ierr);
1767*47c6ae99SBarry Smith 
1768*47c6ae99SBarry Smith   /*
1769*47c6ae99SBarry Smith     For each node in the grid: we get the neighbors in the local (on processor ordering
1770*47c6ae99SBarry Smith     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1771*47c6ae99SBarry Smith     PETSc ordering.
1772*47c6ae99SBarry Smith   */
1773*47c6ae99SBarry Smith   if (!dd->prealloc_only) {
1774*47c6ae99SBarry Smith     ierr = PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);CHKERRQ(ierr);
1775*47c6ae99SBarry Smith     ierr = PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));CHKERRQ(ierr);
1776*47c6ae99SBarry Smith     for (i=xs; i<xs+nx; i++) {
1777*47c6ae99SBarry Smith       istart = DAXPeriodic(wrap) ? -s : (PetscMax(-s,-i));
1778*47c6ae99SBarry Smith       iend   = DAXPeriodic(wrap) ?  s : (PetscMin(s,m-i-1));
1779*47c6ae99SBarry Smith       for (j=ys; j<ys+ny; j++) {
1780*47c6ae99SBarry Smith 	jstart = DAYPeriodic(wrap) ? -s : (PetscMax(-s,-j));
1781*47c6ae99SBarry Smith 	jend   = DAYPeriodic(wrap) ?  s : (PetscMin(s,n-j-1));
1782*47c6ae99SBarry Smith 	for (k=zs; k<zs+nz; k++) {
1783*47c6ae99SBarry Smith 	  kstart = DAZPeriodic(wrap) ? -s : (PetscMax(-s,-k));
1784*47c6ae99SBarry Smith 	  kend   = DAZPeriodic(wrap) ?  s : (PetscMin(s,p-k-1));
1785*47c6ae99SBarry Smith 
1786*47c6ae99SBarry Smith 	  slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1787*47c6ae99SBarry Smith 
1788*47c6ae99SBarry Smith 	  for (l=0; l<nc; l++) {
1789*47c6ae99SBarry Smith 	    cnt  = 0;
1790*47c6ae99SBarry Smith 	    for (ii=istart; ii<iend+1; ii++) {
1791*47c6ae99SBarry Smith 	      for (jj=jstart; jj<jend+1; jj++) {
1792*47c6ae99SBarry Smith 		for (kk=kstart; kk<kend+1; kk++) {
1793*47c6ae99SBarry Smith 		  if (ii || jj || kk) {
1794*47c6ae99SBarry Smith 		    if ((st == DA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1795*47c6ae99SBarry Smith 		      for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
1796*47c6ae99SBarry Smith 			cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1797*47c6ae99SBarry Smith 		    }
1798*47c6ae99SBarry Smith 		  } else {
1799*47c6ae99SBarry Smith 		    if (dfill) {
1800*47c6ae99SBarry Smith 		      for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
1801*47c6ae99SBarry Smith 			cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1802*47c6ae99SBarry Smith 		    } else {
1803*47c6ae99SBarry Smith 		      for (ifill_col=0; ifill_col<nc; ifill_col++)
1804*47c6ae99SBarry Smith 			cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1805*47c6ae99SBarry Smith 		    }
1806*47c6ae99SBarry Smith 		  }
1807*47c6ae99SBarry Smith 		}
1808*47c6ae99SBarry Smith 	      }
1809*47c6ae99SBarry Smith 	    }
1810*47c6ae99SBarry Smith 	    row  = l + nc*(slot);
1811*47c6ae99SBarry Smith 	    ierr = MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);CHKERRQ(ierr);
1812*47c6ae99SBarry Smith 	  }
1813*47c6ae99SBarry Smith 	}
1814*47c6ae99SBarry Smith       }
1815*47c6ae99SBarry Smith     }
1816*47c6ae99SBarry Smith     ierr = PetscFree(values);CHKERRQ(ierr);
1817*47c6ae99SBarry Smith     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1818*47c6ae99SBarry Smith     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1819*47c6ae99SBarry Smith   }
1820*47c6ae99SBarry Smith   ierr = PetscFree(cols);CHKERRQ(ierr);
1821*47c6ae99SBarry Smith   PetscFunctionReturn(0);
1822*47c6ae99SBarry Smith }
1823