Actual source code: fdda.c

  1: 
  2: #include <private/daimpl.h> /*I      "petscdmda.h"     I*/
  3: #include <petscmat.h>         /*I      "petscmat.h"    I*/
  4: #include <private/matimpl.h>

  6: extern PetscErrorCode DMGetColoring_DA_1d_MPIAIJ(DM,ISColoringType,ISColoring *);
  7: extern PetscErrorCode DMGetColoring_DA_2d_MPIAIJ(DM,ISColoringType,ISColoring *);
  8: extern PetscErrorCode DMGetColoring_DA_2d_5pt_MPIAIJ(DM,ISColoringType,ISColoring *);
  9: extern PetscErrorCode DMGetColoring_DA_3d_MPIAIJ(DM,ISColoringType,ISColoring *);

 11: /*
 12:    For ghost i that may be negative or greater than the upper bound this
 13:   maps it into the 0:m-1 range using periodicity
 14: */
 15: #define SetInRange(i,m) ((i < 0) ? m+i:((i >= m) ? i-m:i))

 19: static PetscErrorCode DMDASetBlockFills_Private(PetscInt *dfill,PetscInt w,PetscInt **rfill)
 20: {
 22:   PetscInt       i,j,nz,*fill;

 25:   if (!dfill) return(0);

 27:   /* count number nonzeros */
 28:   nz = 0;
 29:   for (i=0; i<w; i++) {
 30:     for (j=0; j<w; j++) {
 31:       if (dfill[w*i+j]) nz++;
 32:     }
 33:   }
 34:   PetscMalloc((nz + w + 1)*sizeof(PetscInt),&fill);
 35:   /* construct modified CSR storage of nonzero structure */
 36:   nz = w + 1;
 37:   for (i=0; i<w; i++) {
 38:     fill[i] = nz;
 39:     for (j=0; j<w; j++) {
 40:       if (dfill[w*i+j]) {
 41:         fill[nz] = j;
 42:         nz++;
 43:       }
 44:     }
 45:   }
 46:   fill[w] = nz;
 47: 
 48:   *rfill = fill;
 49:   return(0);
 50: }

 54: /*@
 55:     DMDASetBlockFills - Sets the fill pattern in each block for a multi-component problem
 56:     of the matrix returned by DMGetMatrix().

 58:     Logically Collective on DMDA

 60:     Input Parameter:
 61: +   da - the distributed array
 62: .   dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block)
 63: -   ofill - the fill pattern in the off-diagonal blocks


 66:     Level: developer

 68:     Notes: This only makes sense when you are doing multicomponent problems but using the
 69:        MPIAIJ matrix format

 71:            The format for dfill and ofill is a 2 dimensional dof by dof matrix with 1 entries
 72:        representing coupling and 0 entries for missing coupling. For example 
 73: $             dfill[9] = {1, 0, 0,
 74: $                         1, 1, 0,
 75: $                         0, 1, 1} 
 76:        means that row 0 is coupled with only itself in the diagonal block, row 1 is coupled with 
 77:        itself and row 0 (in the diagonal block) and row 2 is coupled with itself and row 1 (in the 
 78:        diagonal block).

 80:      DMDASetGetMatrix() allows you to provide general code for those more complicated nonzero patterns then
 81:      can be represented in the dfill, ofill format

 83:    Contributed by Glenn Hammond

 85: .seealso DMGetMatrix(), DMDASetGetMatrix(), DMDASetMatPreallocateOnly()

 87: @*/
 88: PetscErrorCode  DMDASetBlockFills(DM da,PetscInt *dfill,PetscInt *ofill)
 89: {
 90:   DM_DA          *dd = (DM_DA*)da->data;

 94:   DMDASetBlockFills_Private(dfill,dd->w,&dd->dfill);
 95:   DMDASetBlockFills_Private(ofill,dd->w,&dd->ofill);
 96:   return(0);
 97: }


102: PetscErrorCode  DMGetColoring_DA(DM da,ISColoringType ctype,const MatType mtype,ISColoring *coloring)
103: {
104:   PetscErrorCode   ierr;
105:   PetscInt         dim,m,n,p,nc;
106:   DMDABoundaryType bx,by,bz;
107:   MPI_Comm         comm;
108:   PetscMPIInt      size;
109:   PetscBool        isBAIJ;
110:   DM_DA            *dd = (DM_DA*)da->data;

113:   /*
114:                                   m
115:           ------------------------------------------------------
116:          |                                                     |
117:          |                                                     |
118:          |               ----------------------                |
119:          |               |                    |                |
120:       n  |           yn  |                    |                |
121:          |               |                    |                |
122:          |               .---------------------                |
123:          |             (xs,ys)     xn                          |
124:          |            .                                        |
125:          |         (gxs,gys)                                   |
126:          |                                                     |
127:           -----------------------------------------------------
128:   */

130:   /*     
131:          nc - number of components per grid point 
132:          col - number of colors needed in one direction for single component problem
133:   
134:   */
135:   DMDAGetInfo(da,&dim,0,0,0,&m,&n,&p,&nc,0,&bx,&by,&bz,0);

137:   PetscObjectGetComm((PetscObject)da,&comm);
138:   MPI_Comm_size(comm,&size);
139:   if (ctype == IS_COLORING_GHOSTED){
140:     if (size == 1) {
141:       ctype = IS_COLORING_GLOBAL;
142:     } else if (dim > 1){
143:       if ((m==1 && bx == DMDA_BOUNDARY_PERIODIC) || (n==1 && by == DMDA_BOUNDARY_PERIODIC) || (p==1 && bz == DMDA_BOUNDARY_PERIODIC)){
144:         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");
145:       }
146:     }
147:   }

149:   /* Tell the DMDA it has 1 degree of freedom per grid point so that the coloring for BAIJ 
150:      matrices is for the blocks, not the individual matrix elements  */
151:   PetscStrcmp(mtype,MATBAIJ,&isBAIJ);
152:   if (!isBAIJ) {PetscStrcmp(mtype,MATMPIBAIJ,&isBAIJ);}
153:   if (!isBAIJ) {PetscStrcmp(mtype,MATSEQBAIJ,&isBAIJ);}
154:   if (isBAIJ) {
155:     dd->w = 1;
156:     dd->xs = dd->xs/nc;
157:     dd->xe = dd->xe/nc;
158:     dd->Xs = dd->Xs/nc;
159:     dd->Xe = dd->Xe/nc;
160:   }

162:   /*
163:      We do not provide a getcoloring function in the DMDA operations because 
164:    the basic DMDA does not know about matrices. We think of DMDA as being more 
165:    more low-level then matrices.
166:   */
167:   if (dim == 1) {
168:     DMGetColoring_DA_1d_MPIAIJ(da,ctype,coloring);
169:   } else if (dim == 2) {
170:      DMGetColoring_DA_2d_MPIAIJ(da,ctype,coloring);
171:   } else if (dim == 3) {
172:      DMGetColoring_DA_3d_MPIAIJ(da,ctype,coloring);
173:   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not done for %D dimension, send us mail petsc-maint@mcs.anl.gov for code",dim);
174:   if (isBAIJ) {
175:     dd->w = nc;
176:     dd->xs = dd->xs*nc;
177:     dd->xe = dd->xe*nc;
178:     dd->Xs = dd->Xs*nc;
179:     dd->Xe = dd->Xe*nc;
180:   }
181:   return(0);
182: }

184: /* ---------------------------------------------------------------------------------*/

188: PetscErrorCode DMGetColoring_DA_2d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
189: {
190:   PetscErrorCode         ierr;
191:   PetscInt               xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,M,N,dim,s,k,nc,col;
192:   PetscInt               ncolors;
193:   MPI_Comm               comm;
194:   DMDABoundaryType       bx,by;
195:   DMDAStencilType        st;
196:   ISColoringValue        *colors;
197:   DM_DA                  *dd = (DM_DA*)da->data;

200:   /*     
201:          nc - number of components per grid point 
202:          col - number of colors needed in one direction for single component problem
203:   
204:   */
205:   DMDAGetInfo(da,&dim,&m,&n,0,&M,&N,0,&nc,&s,&bx,&by,0,&st);
206:   col    = 2*s + 1;
207:   DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);
208:   DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);
209:   PetscObjectGetComm((PetscObject)da,&comm);

211:   /* special case as taught to us by Paul Hovland */
212:   if (st == DMDA_STENCIL_STAR && s == 1) {
213:     DMGetColoring_DA_2d_5pt_MPIAIJ(da,ctype,coloring);
214:   } else {

216:     if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){
217:       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X (%d) is divisible\n\
218:                  by 2*stencil_width + 1 (%d)\n", m, col);
219:     }
220:     if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){
221:       SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y (%d) is divisible\n\
222:                  by 2*stencil_width + 1 (%d)\n", n, col);
223:     }
224:     if (ctype == IS_COLORING_GLOBAL) {
225:       if (!dd->localcoloring) {
226:         PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);
227:         ii = 0;
228:         for (j=ys; j<ys+ny; j++) {
229:           for (i=xs; i<xs+nx; i++) {
230:             for (k=0; k<nc; k++) {
231:               colors[ii++] = k + nc*((i % col) + col*(j % col));
232:             }
233:           }
234:         }
235:         ncolors = nc + nc*(col-1 + col*(col-1));
236:         ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);
237:       }
238:       *coloring = dd->localcoloring;
239:     } else if (ctype == IS_COLORING_GHOSTED) {
240:       if (!dd->ghostedcoloring) {
241:         PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);
242:         ii = 0;
243:         for (j=gys; j<gys+gny; j++) {
244:           for (i=gxs; i<gxs+gnx; i++) {
245:             for (k=0; k<nc; k++) {
246:               /* the complicated stuff is to handle periodic boundaries */
247:               colors[ii++] = k + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col));
248:             }
249:           }
250:         }
251:         ncolors = nc + nc*(col - 1 + col*(col-1));
252:         ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);
253:         /* PetscIntView(ncolors,(PetscInt *)colors,0); */

255:         ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);
256:       }
257:       *coloring = dd->ghostedcoloring;
258:     } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
259:   }
260:   ISColoringReference(*coloring);
261:   return(0);
262: }

264: /* ---------------------------------------------------------------------------------*/

268: PetscErrorCode DMGetColoring_DA_3d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
269: {
270:   PetscErrorCode    ierr;
271:   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;
272:   PetscInt          ncolors;
273:   MPI_Comm          comm;
274:   DMDABoundaryType  bx,by,bz;
275:   DMDAStencilType   st;
276:   ISColoringValue   *colors;
277:   DM_DA             *dd = (DM_DA*)da->data;

280:   /*     
281:          nc - number of components per grid point 
282:          col - number of colors needed in one direction for single component problem
283:   
284:   */
285:   DMDAGetInfo(da,&dim,&m,&n,&p,&M,&N,&P,&nc,&s,&bx,&by,&bz,&st);
286:   col    = 2*s + 1;
287:   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){
288:     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
289:                  by 2*stencil_width + 1\n");
290:   }
291:   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){
292:     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
293:                  by 2*stencil_width + 1\n");
294:   }
295:   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)){
296:     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
297:                  by 2*stencil_width + 1\n");
298:   }

300:   DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);
301:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);
302:   PetscObjectGetComm((PetscObject)da,&comm);

304:   /* create the coloring */
305:   if (ctype == IS_COLORING_GLOBAL) {
306:     if (!dd->localcoloring) {
307:       PetscMalloc(nc*nx*ny*nz*sizeof(ISColoringValue),&colors);
308:       ii = 0;
309:       for (k=zs; k<zs+nz; k++) {
310:         for (j=ys; j<ys+ny; j++) {
311:           for (i=xs; i<xs+nx; i++) {
312:             for (l=0; l<nc; l++) {
313:               colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
314:             }
315:           }
316:         }
317:       }
318:       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
319:       ISColoringCreate(comm,ncolors,nc*nx*ny*nz,colors,&dd->localcoloring);
320:     }
321:     *coloring = dd->localcoloring;
322:   } else if (ctype == IS_COLORING_GHOSTED) {
323:     if (!dd->ghostedcoloring) {
324:       PetscMalloc(nc*gnx*gny*gnz*sizeof(ISColoringValue),&colors);
325:       ii = 0;
326:       for (k=gzs; k<gzs+gnz; k++) {
327:         for (j=gys; j<gys+gny; j++) {
328:           for (i=gxs; i<gxs+gnx; i++) {
329:             for (l=0; l<nc; l++) {
330:               /* the complicated stuff is to handle periodic boundaries */
331:               colors[ii++] = l + nc*((SetInRange(i,m) % col) + col*(SetInRange(j,n) % col) + col*col*(SetInRange(k,p) % col));
332:             }
333:           }
334:         }
335:       }
336:       ncolors = nc + nc*(col-1 + col*(col-1)+ col*col*(col-1));
337:       ISColoringCreate(comm,ncolors,nc*gnx*gny*gnz,colors,&dd->ghostedcoloring);
338:       ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);
339:     }
340:     *coloring = dd->ghostedcoloring;
341:   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
342:   ISColoringReference(*coloring);
343:   return(0);
344: }

346: /* ---------------------------------------------------------------------------------*/

350: PetscErrorCode DMGetColoring_DA_1d_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
351: {
352:   PetscErrorCode    ierr;
353:   PetscInt          xs,nx,i,i1,gxs,gnx,l,m,M,dim,s,nc,col;
354:   PetscInt          ncolors;
355:   MPI_Comm          comm;
356:   DMDABoundaryType  bx;
357:   ISColoringValue   *colors;
358:   DM_DA             *dd = (DM_DA*)da->data;

361:   /*     
362:          nc - number of components per grid point 
363:          col - number of colors needed in one direction for single component problem
364:   
365:   */
366:   DMDAGetInfo(da,&dim,&m,0,0,&M,0,0,&nc,&s,&bx,0,0,0);
367:   col    = 2*s + 1;

369:   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)) {
370:     SETERRQ2(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points %d is divisible\n\
371:                  by 2*stencil_width + 1 %d\n",(int)m,(int)col);
372:   }

374:   DMDAGetCorners(da,&xs,0,0,&nx,0,0);
375:   DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);
376:   PetscObjectGetComm((PetscObject)da,&comm);

378:   /* create the coloring */
379:   if (ctype == IS_COLORING_GLOBAL) {
380:     if (!dd->localcoloring) {
381:       PetscMalloc(nc*nx*sizeof(ISColoringValue),&colors);
382:       i1 = 0;
383:       for (i=xs; i<xs+nx; i++) {
384:         for (l=0; l<nc; l++) {
385:           colors[i1++] = l + nc*(i % col);
386:         }
387:       }
388:       ncolors = nc + nc*(col-1);
389:       ISColoringCreate(comm,ncolors,nc*nx,colors,&dd->localcoloring);
390:     }
391:     *coloring = dd->localcoloring;
392:   } else if (ctype == IS_COLORING_GHOSTED) {
393:     if (!dd->ghostedcoloring) {
394:       PetscMalloc(nc*gnx*sizeof(ISColoringValue),&colors);
395:       i1 = 0;
396:       for (i=gxs; i<gxs+gnx; i++) {
397:         for (l=0; l<nc; l++) {
398:           /* the complicated stuff is to handle periodic boundaries */
399:           colors[i1++] = l + nc*(SetInRange(i,m) % col);
400:         }
401:       }
402:       ncolors = nc + nc*(col-1);
403:       ISColoringCreate(comm,ncolors,nc*gnx,colors,&dd->ghostedcoloring);
404:       ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);
405:     }
406:     *coloring = dd->ghostedcoloring;
407:   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
408:   ISColoringReference(*coloring);
409:   return(0);
410: }

414: PetscErrorCode DMGetColoring_DA_2d_5pt_MPIAIJ(DM da,ISColoringType ctype,ISColoring *coloring)
415: {
416:   PetscErrorCode    ierr;
417:   PetscInt          xs,ys,nx,ny,i,j,ii,gxs,gys,gnx,gny,m,n,dim,s,k,nc;
418:   PetscInt          ncolors;
419:   MPI_Comm          comm;
420:   DMDABoundaryType  bx,by;
421:   ISColoringValue   *colors;
422:   DM_DA             *dd = (DM_DA*)da->data;

425:   /*     
426:          nc - number of components per grid point 
427:          col - number of colors needed in one direction for single component problem
428:   
429:   */
430:   DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,0);
431:   DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);
432:   DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);
433:   PetscObjectGetComm((PetscObject)da,&comm);

435:   if (bx == DMDA_BOUNDARY_PERIODIC && (m % 5)){
436:     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
437:                  by 5\n");
438:   }
439:   if (by == DMDA_BOUNDARY_PERIODIC && (n % 5)){
440:     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
441:                  by 5\n");
442:   }

444:   /* create the coloring */
445:   if (ctype == IS_COLORING_GLOBAL) {
446:     if (!dd->localcoloring) {
447:       PetscMalloc(nc*nx*ny*sizeof(ISColoringValue),&colors);
448:       ii = 0;
449:       for (j=ys; j<ys+ny; j++) {
450:         for (i=xs; i<xs+nx; i++) {
451:           for (k=0; k<nc; k++) {
452:             colors[ii++] = k + nc*((3*j+i) % 5);
453:           }
454:         }
455:       }
456:       ncolors = 5*nc;
457:       ISColoringCreate(comm,ncolors,nc*nx*ny,colors,&dd->localcoloring);
458:     }
459:     *coloring = dd->localcoloring;
460:   } else if (ctype == IS_COLORING_GHOSTED) {
461:     if (!dd->ghostedcoloring) {
462:       PetscMalloc(nc*gnx*gny*sizeof(ISColoringValue),&colors);
463:       ii = 0;
464:       for (j=gys; j<gys+gny; j++) {
465:         for (i=gxs; i<gxs+gnx; i++) {
466:           for (k=0; k<nc; k++) {
467:             colors[ii++] = k + nc*((3*SetInRange(j,n) + SetInRange(i,m)) % 5);
468:           }
469:         }
470:       }
471:       ncolors = 5*nc;
472:       ISColoringCreate(comm,ncolors,nc*gnx*gny,colors,&dd->ghostedcoloring);
473:       ISColoringSetType(dd->ghostedcoloring,IS_COLORING_GHOSTED);
474:     }
475:     *coloring = dd->ghostedcoloring;
476:   } else SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Unknown ISColoringType %d",(int)ctype);
477:   return(0);
478: }

480: /* =========================================================================== */
481: extern PetscErrorCode DMGetMatrix_DA_1d_MPIAIJ(DM,Mat);
482: extern PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ(DM,Mat);
483: extern PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ_Fill(DM,Mat);
484: extern PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ(DM,Mat);
485: extern PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ_Fill(DM,Mat);
486: extern PetscErrorCode DMGetMatrix_DA_2d_MPIBAIJ(DM,Mat);
487: extern PetscErrorCode DMGetMatrix_DA_3d_MPIBAIJ(DM,Mat);
488: extern PetscErrorCode DMGetMatrix_DA_2d_MPISBAIJ(DM,Mat);
489: extern PetscErrorCode DMGetMatrix_DA_3d_MPISBAIJ(DM,Mat);

493: /*@
494:    MatSetDM - Sets the DMDA that is to be used by the HYPRE_StructMatrix PETSc matrix

496:    Logically Collective on Mat

498:    Input Parameters:
499: +  mat - the matrix
500: -  da - the da

502:    Level: intermediate

504: @*/
505: PetscErrorCode  MatSetDM(Mat mat,DM da)
506: {

512:   PetscTryMethod(mat,"MatSetDM_C",(Mat,DM),(mat,da));
513:   return(0);
514: }

516: EXTERN_C_BEGIN
519: PetscErrorCode  MatView_MPI_DA(Mat A,PetscViewer viewer)
520: {
521:   DM             da;
523:   const char     *prefix;
524:   Mat            Anatural;
525:   AO             ao;
526:   PetscInt       rstart,rend,*petsc,i;
527:   IS             is;
528:   MPI_Comm       comm;
529:   PetscViewerFormat format;

532:   /* Check whether we are just printing info, in which case MatView() already viewed everything we wanted to view */
533:   PetscViewerGetFormat(viewer,&format);
534:   if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) return(0);

536:   PetscObjectGetComm((PetscObject)A,&comm);
537:   PetscObjectQuery((PetscObject)A,"DM",(PetscObject*)&da);
538:   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");

540:   DMDAGetAO(da,&ao);
541:   MatGetOwnershipRange(A,&rstart,&rend);
542:   PetscMalloc((rend-rstart)*sizeof(PetscInt),&petsc);
543:   for (i=rstart; i<rend; i++) petsc[i-rstart] = i;
544:   AOApplicationToPetsc(ao,rend-rstart,petsc);
545:   ISCreateGeneral(comm,rend-rstart,petsc,PETSC_OWN_POINTER,&is);

547:   /* call viewer on natural ordering */
548:   MatGetSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&Anatural);
549:   ISDestroy(&is);
550:   PetscObjectGetOptionsPrefix((PetscObject)A,&prefix);
551:   PetscObjectSetOptionsPrefix((PetscObject)Anatural,prefix);
552:   PetscObjectSetName((PetscObject)Anatural,((PetscObject)A)->name);
553:   MatView(Anatural,viewer);
554:   MatDestroy(&Anatural);
555:   return(0);
556: }
557: EXTERN_C_END

559: EXTERN_C_BEGIN
562: PetscErrorCode  MatLoad_MPI_DA(Mat A,PetscViewer viewer)
563: {
564:   DM             da;
566:   Mat            Anatural,Aapp;
567:   AO             ao;
568:   PetscInt       rstart,rend,*app,i;
569:   IS             is;
570:   MPI_Comm       comm;

573:   PetscObjectGetComm((PetscObject)A,&comm);
574:   PetscObjectQuery((PetscObject)A,"DM",(PetscObject*)&da);
575:   if (!da) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"Matrix not generated from a DMDA");

577:   /* Load the matrix in natural ordering */
578:   MatCreate(((PetscObject)A)->comm,&Anatural);
579:   MatSetType(Anatural,((PetscObject)A)->type_name);
580:   MatSetSizes(Anatural,A->rmap->n,A->cmap->n,A->rmap->N,A->cmap->N);
581:   MatLoad(Anatural,viewer);

583:   /* Map natural ordering to application ordering and create IS */
584:   DMDAGetAO(da,&ao);
585:   MatGetOwnershipRange(Anatural,&rstart,&rend);
586:   PetscMalloc((rend-rstart)*sizeof(PetscInt),&app);
587:   for (i=rstart; i<rend; i++) app[i-rstart] = i;
588:   AOPetscToApplication(ao,rend-rstart,app);
589:   ISCreateGeneral(comm,rend-rstart,app,PETSC_OWN_POINTER,&is);

591:   /* Do permutation and replace header */
592:   MatGetSubMatrix(Anatural,is,is,MAT_INITIAL_MATRIX,&Aapp);
593:   MatHeaderReplace(A,Aapp);
594:   ISDestroy(&is);
595:   MatDestroy(&Anatural);
596:   return(0);
597: }
598: EXTERN_C_END

602: PetscErrorCode  DMGetMatrix_DA(DM da, const MatType mtype,Mat *J)
603: {
605:   PetscInt       dim,dof,nx,ny,nz,dims[3],starts[3],M,N,P;
606:   Mat            A;
607:   MPI_Comm       comm;
608:   const MatType  Atype;
609:   void           (*aij)(void)=PETSC_NULL,(*baij)(void)=PETSC_NULL,(*sbaij)(void)=PETSC_NULL;
610:   MatType        ttype[256];
611:   PetscBool      flg;
612:   PetscMPIInt    size;
613:   DM_DA          *dd = (DM_DA*)da->data;

616: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
617:   MatInitializePackage(PETSC_NULL);
618: #endif
619:   if (!mtype) mtype = MATAIJ;
620:   PetscStrcpy((char*)ttype,mtype);
621:   PetscOptionsBegin(((PetscObject)da)->comm,((PetscObject)da)->prefix,"DMDA options","Mat");
622:   PetscOptionsList("-da_mat_type","Matrix type","MatSetType",MatList,mtype,(char*)ttype,256,&flg);
623:   PetscOptionsEnd();

625:   /*
626:                                   m
627:           ------------------------------------------------------
628:          |                                                     |
629:          |                                                     |
630:          |               ----------------------                |
631:          |               |                    |                |
632:       n  |           ny  |                    |                |
633:          |               |                    |                |
634:          |               .---------------------                |
635:          |             (xs,ys)     nx                          |
636:          |            .                                        |
637:          |         (gxs,gys)                                   |
638:          |                                                     |
639:           -----------------------------------------------------
640:   */

642:   /*     
643:          nc - number of components per grid point 
644:          col - number of colors needed in one direction for single component problem
645:   
646:   */
647:   DMDAGetInfo(da,&dim,&M,&N,&P,0,0,0,&dof,0,0,0,0,0);
648:   DMDAGetCorners(da,0,0,0,&nx,&ny,&nz);
649:   PetscObjectGetComm((PetscObject)da,&comm);
650:   MatCreate(comm,&A);
651:   MatSetSizes(A,dof*nx*ny*nz,dof*nx*ny*nz,dof*M*N*P,dof*M*N*P);
652:   MatSetType(A,(const MatType)ttype);
653:   MatSetDM(A,da);
654:   MatSetFromOptions(A);
655:   MatGetType(A,&Atype);
656:   /*
657:      We do not provide a getmatrix function in the DMDA operations because 
658:    the basic DMDA does not know about matrices. We think of DMDA as being more 
659:    more low-level than matrices. This is kind of cheating but, cause sometimes 
660:    we think of DMDA has higher level than matrices.

662:      We could switch based on Atype (or mtype), but we do not since the
663:    specialized setting routines depend only the particular preallocation
664:    details of the matrix, not the type itself.
665:   */
666:   PetscObjectQueryFunction((PetscObject)A,"MatMPIAIJSetPreallocation_C",&aij);
667:   if (!aij) {
668:     PetscObjectQueryFunction((PetscObject)A,"MatSeqAIJSetPreallocation_C",&aij);
669:   }
670:   if (!aij) {
671:     PetscObjectQueryFunction((PetscObject)A,"MatMPIBAIJSetPreallocation_C",&baij);
672:     if (!baij) {
673:       PetscObjectQueryFunction((PetscObject)A,"MatSeqBAIJSetPreallocation_C",&baij);
674:     }
675:     if (!baij){
676:       PetscObjectQueryFunction((PetscObject)A,"MatMPISBAIJSetPreallocation_C",&sbaij);
677:       if (!sbaij) {
678:         PetscObjectQueryFunction((PetscObject)A,"MatSeqSBAIJSetPreallocation_C",&sbaij);
679:       }
680:     }
681:   }
682:   if (aij) {
683:     if (dim == 1) {
684:       DMGetMatrix_DA_1d_MPIAIJ(da,A);
685:     } else if (dim == 2) {
686:       if (dd->ofill) {
687:         DMGetMatrix_DA_2d_MPIAIJ_Fill(da,A);
688:       } else {
689:         DMGetMatrix_DA_2d_MPIAIJ(da,A);
690:       }
691:     } else if (dim == 3) {
692:       if (dd->ofill) {
693:         DMGetMatrix_DA_3d_MPIAIJ_Fill(da,A);
694:       } else {
695:         DMGetMatrix_DA_3d_MPIAIJ(da,A);
696:       }
697:     }
698:   } else if (baij) {
699:     if (dim == 2) {
700:       DMGetMatrix_DA_2d_MPIBAIJ(da,A);
701:     } else if (dim == 3) {
702:       DMGetMatrix_DA_3d_MPIBAIJ(da,A);
703:     } else {
704:       SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \
705:                "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim);
706:     }
707:   } else if (sbaij) {
708:     if (dim == 2) {
709:       DMGetMatrix_DA_2d_MPISBAIJ(da,A);
710:     } else if (dim == 3) {
711:       DMGetMatrix_DA_3d_MPISBAIJ(da,A);
712:     } else {
713:       SETERRQ3(((PetscObject)da)->comm,PETSC_ERR_SUP,"Not implemented for %D dimension and Matrix Type: %s in %D dimension!\n" \
714:                "Send mail to petsc-maint@mcs.anl.gov for code",dim,Atype,dim);
715:     }
716:   } else {
717:     ISLocalToGlobalMapping ltog,ltogb;
718:     DMGetLocalToGlobalMapping(da,&ltog);
719:     DMGetLocalToGlobalMappingBlock(da,&ltogb);
720:     MatSetLocalToGlobalMapping(A,ltog,ltog);
721:     MatSetLocalToGlobalMappingBlock(A,ltogb,ltogb);
722:   }
723:   DMDAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);
724:   MatSetStencil(A,dim,dims,starts,dof);
725:   PetscObjectCompose((PetscObject)A,"DM",(PetscObject)da);
726:   MPI_Comm_size(comm,&size);
727:   if (size > 1) {
728:     /* change viewer to display matrix in natural ordering */
729:     MatShellSetOperation(A, MATOP_VIEW, (void (*)(void)) MatView_MPI_DA);
730:     MatShellSetOperation(A, MATOP_LOAD, (void (*)(void)) MatLoad_MPI_DA);
731:   }
732:   *J = A;
733:   return(0);
734: }

736: /* ---------------------------------------------------------------------------------*/
739: PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ(DM da,Mat J)
740: {
741:   PetscErrorCode         ierr;
742:   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;
743:   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
744:   MPI_Comm               comm;
745:   PetscScalar            *values;
746:   DMDABoundaryType       bx,by;
747:   ISLocalToGlobalMapping ltog,ltogb;
748:   DMDAStencilType        st;

751:   /*     
752:          nc - number of components per grid point 
753:          col - number of colors needed in one direction for single component problem
754:   
755:   */
756:   DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);
757:   col = 2*s + 1;
758:   DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);
759:   DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);
760:   PetscObjectGetComm((PetscObject)da,&comm);

762:   PetscMalloc2(nc,PetscInt,&rows,col*col*nc*nc,PetscInt,&cols);
763:   DMGetLocalToGlobalMapping(da,&ltog);
764:   DMGetLocalToGlobalMappingBlock(da,&ltogb);
765: 
766:   /* determine the matrix preallocation information */
767:   MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);
768:   for (i=xs; i<xs+nx; i++) {

770:     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
771:     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));

773:     for (j=ys; j<ys+ny; j++) {
774:       slot = i - gxs + gnx*(j - gys);

776:       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
777:       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));

779:       cnt  = 0;
780:       for (k=0; k<nc; k++) {
781:         for (l=lstart; l<lend+1; l++) {
782:           for (p=pstart; p<pend+1; p++) {
783:             if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
784:               cols[cnt++]  = k + nc*(slot + gnx*l + p);
785:             }
786:           }
787:         }
788:         rows[k] = k + nc*(slot);
789:       }
790:       MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);
791:     }
792:   }
793:   MatSeqAIJSetPreallocation(J,0,dnz);
794:   MatMPIAIJSetPreallocation(J,0,dnz,0,onz);
795:   MatSetBlockSize(J,nc);
796:   MatPreallocateFinalize(dnz,onz);

798:   MatSetLocalToGlobalMapping(J,ltog,ltog);
799:   MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);

801:   /*
802:     For each node in the grid: we get the neighbors in the local (on processor ordering
803:     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
804:     PETSc ordering.
805:   */
806:   if (!da->prealloc_only) {
807:     PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);
808:     PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));
809:     for (i=xs; i<xs+nx; i++) {
810: 
811:       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
812:       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
813: 
814:       for (j=ys; j<ys+ny; j++) {
815:         slot = i - gxs + gnx*(j - gys);
816: 
817:         lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
818:         lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));

820:         cnt  = 0;
821:         for (k=0; k<nc; k++) {
822:           for (l=lstart; l<lend+1; l++) {
823:             for (p=pstart; p<pend+1; p++) {
824:               if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star have either l = 0 or p = 0 */
825:                 cols[cnt++]  = k + nc*(slot + gnx*l + p);
826:               }
827:             }
828:           }
829:           rows[k]      = k + nc*(slot);
830:         }
831:         MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);
832:       }
833:     }
834:     PetscFree(values);
835:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
836:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
837:   }
838:   PetscFree2(rows,cols);
839:   return(0);
840: }

844: PetscErrorCode DMGetMatrix_DA_2d_MPIAIJ_Fill(DM da,Mat J)
845: {
846:   PetscErrorCode         ierr;
847:   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
848:   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p;
849:   PetscInt               lstart,lend,pstart,pend,*dnz,*onz;
850:   DM_DA                  *dd = (DM_DA*)da->data;
851:   PetscInt               ifill_col,*ofill = dd->ofill, *dfill = dd->dfill;
852:   MPI_Comm               comm;
853:   PetscScalar            *values;
854:   DMDABoundaryType       bx,by;
855:   ISLocalToGlobalMapping ltog,ltogb;
856:   DMDAStencilType        st;

859:   /*     
860:          nc - number of components per grid point 
861:          col - number of colors needed in one direction for single component problem
862:   
863:   */
864:   DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);
865:   col = 2*s + 1;
866:   DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);
867:   DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);
868:   PetscObjectGetComm((PetscObject)da,&comm);

870:   PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);
871:   DMGetLocalToGlobalMapping(da,&ltog);
872:   DMGetLocalToGlobalMappingBlock(da,&ltogb);
873: 
874:   /* determine the matrix preallocation information */
875:   MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);
876:   for (i=xs; i<xs+nx; i++) {

878:     pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
879:     pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));

881:     for (j=ys; j<ys+ny; j++) {
882:       slot = i - gxs + gnx*(j - gys);

884:       lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
885:       lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));

887:       for (k=0; k<nc; k++) {
888:         cnt  = 0;
889:         for (l=lstart; l<lend+1; l++) {
890:           for (p=pstart; p<pend+1; p++) {
891:             if (l || p) {
892:               if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
893:                 for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
894:                   cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
895:               }
896:             } else {
897:               if (dfill) {
898:                 for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
899:                   cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
900:               } else {
901:                 for (ifill_col=0; ifill_col<nc; ifill_col++)
902:                   cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
903:               }
904:             }
905:           }
906:         }
907:         row = k + nc*(slot);
908:         MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);
909:       }
910:     }
911:   }
912:   MatSeqAIJSetPreallocation(J,0,dnz);
913:   MatMPIAIJSetPreallocation(J,0,dnz,0,onz);
914:   MatPreallocateFinalize(dnz,onz);
915:   MatSetLocalToGlobalMapping(J,ltog,ltog);
916:   MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);

918:   /*
919:     For each node in the grid: we get the neighbors in the local (on processor ordering
920:     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
921:     PETSc ordering.
922:   */
923:   if (!da->prealloc_only) {
924:     PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);
925:     PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));
926:     for (i=xs; i<xs+nx; i++) {
927: 
928:       pstart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
929:       pend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
930: 
931:       for (j=ys; j<ys+ny; j++) {
932:         slot = i - gxs + gnx*(j - gys);
933: 
934:         lstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
935:         lend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));

937:         for (k=0; k<nc; k++) {
938:           cnt  = 0;
939:           for (l=lstart; l<lend+1; l++) {
940:             for (p=pstart; p<pend+1; p++) {
941:               if (l || p) {
942:                 if ((st == DMDA_STENCIL_BOX) || (!l || !p)) {  /* entries on star */
943:                   for (ifill_col=ofill[k]; ifill_col<ofill[k+1]; ifill_col++)
944:                     cols[cnt++]  = ofill[ifill_col] + nc*(slot + gnx*l + p);
945:                 }
946:               } else {
947:                 if (dfill) {
948:                   for (ifill_col=dfill[k]; ifill_col<dfill[k+1]; ifill_col++)
949:                     cols[cnt++]  = dfill[ifill_col] + nc*(slot + gnx*l + p);
950:                 } else {
951:                   for (ifill_col=0; ifill_col<nc; ifill_col++)
952:                     cols[cnt++]  = ifill_col + nc*(slot + gnx*l + p);
953:                 }
954:               }
955:             }
956:           }
957:           row  = k + nc*(slot);
958:           MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);
959:         }
960:       }
961:     }
962:     PetscFree(values);
963:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
964:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
965:   }
966:   PetscFree(cols);
967:   return(0);
968: }

970: /* ---------------------------------------------------------------------------------*/

974: PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ(DM da,Mat J)
975: {
976:   PetscErrorCode         ierr;
977:   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
978:   PetscInt               m,n,dim,s,*cols = PETSC_NULL,k,nc,*rows = PETSC_NULL,col,cnt,l,p,*dnz = PETSC_NULL,*onz = PETSC_NULL;
979:   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
980:   MPI_Comm               comm;
981:   PetscScalar            *values;
982:   DMDABoundaryType       bx,by,bz;
983:   ISLocalToGlobalMapping ltog,ltogb;
984:   DMDAStencilType        st;

987:   /*     
988:          nc - number of components per grid point 
989:          col - number of colors needed in one direction for single component problem
990:   
991:   */
992:   DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);
993:   col    = 2*s + 1;

995:   DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);
996:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);
997:   PetscObjectGetComm((PetscObject)da,&comm);

999:   PetscMalloc2(nc,PetscInt,&rows,col*col*col*nc*nc,PetscInt,&cols);
1000:   DMGetLocalToGlobalMapping(da,&ltog);
1001:   DMGetLocalToGlobalMappingBlock(da,&ltogb);

1003:   /* determine the matrix preallocation information */
1004:   MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);
1005:   for (i=xs; i<xs+nx; i++) {
1006:     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1007:     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1008:     for (j=ys; j<ys+ny; j++) {
1009:       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1010:       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1011:       for (k=zs; k<zs+nz; k++) {
1012:         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1013:         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
1014: 
1015:         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1016: 
1017:         cnt  = 0;
1018:         for (l=0; l<nc; l++) {
1019:           for (ii=istart; ii<iend+1; ii++) {
1020:             for (jj=jstart; jj<jend+1; jj++) {
1021:               for (kk=kstart; kk<kend+1; kk++) {
1022:                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1023:                   cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1024:                 }
1025:               }
1026:             }
1027:           }
1028:           rows[l] = l + nc*(slot);
1029:         }
1030:         MatPreallocateSetLocal(ltog,nc,rows,ltog,cnt,cols,dnz,onz);
1031:       }
1032:     }
1033:   }
1034:   MatSeqAIJSetPreallocation(J,0,dnz);
1035:   MatMPIAIJSetPreallocation(J,0,dnz,0,onz);
1036:   MatPreallocateFinalize(dnz,onz);
1037:   MatSetBlockSize(J,nc);
1038:   MatSetLocalToGlobalMapping(J,ltog,ltog);
1039:   MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);

1041:   /*
1042:     For each node in the grid: we get the neighbors in the local (on processor ordering
1043:     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1044:     PETSc ordering.
1045:   */
1046:   if (!da->prealloc_only) {
1047:     PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);
1048:     PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));
1049:     for (i=xs; i<xs+nx; i++) {
1050:       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1051:       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1052:       for (j=ys; j<ys+ny; j++) {
1053:         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1054:         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1055:         for (k=zs; k<zs+nz; k++) {
1056:           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1057:           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
1058: 
1059:           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1060: 
1061:           cnt  = 0;
1062:           for (l=0; l<nc; l++) {
1063:             for (ii=istart; ii<iend+1; ii++) {
1064:               for (jj=jstart; jj<jend+1; jj++) {
1065:                 for (kk=kstart; kk<kend+1; kk++) {
1066:                   if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1067:                     cols[cnt++]  = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1068:                   }
1069:                 }
1070:               }
1071:             }
1072:             rows[l]      = l + nc*(slot);
1073:           }
1074:           MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);
1075:         }
1076:       }
1077:     }
1078:     PetscFree(values);
1079:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
1080:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
1081:   }
1082:   PetscFree2(rows,cols);
1083:   return(0);
1084: }

1086: /* ---------------------------------------------------------------------------------*/

1090: PetscErrorCode DMGetMatrix_DA_1d_MPIAIJ(DM da,Mat J)
1091: {
1092:   PetscErrorCode         ierr;
1093:   PetscInt               xs,nx,i,i1,slot,gxs,gnx;
1094:   PetscInt               m,dim,s,*cols = PETSC_NULL,nc,*rows = PETSC_NULL,col,cnt,l;
1095:   PetscInt               istart,iend;
1096:   PetscScalar            *values;
1097:   DMDABoundaryType       bx;
1098:   ISLocalToGlobalMapping ltog,ltogb;

1101:   /*     
1102:          nc - number of components per grid point 
1103:          col - number of colors needed in one direction for single component problem
1104:   
1105:   */
1106:   DMDAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&bx,0,0,0);
1107:   col    = 2*s + 1;

1109:   DMDAGetCorners(da,&xs,0,0,&nx,0,0);
1110:   DMDAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);

1112:   MatSeqAIJSetPreallocation(J,col*nc,0);
1113:   MatMPIAIJSetPreallocation(J,col*nc,0,col*nc,0);
1114:   MatSetBlockSize(J,nc);
1115:   PetscMalloc2(nc,PetscInt,&rows,col*nc*nc,PetscInt,&cols);
1116: 
1117:   DMGetLocalToGlobalMapping(da,&ltog);
1118:   DMGetLocalToGlobalMappingBlock(da,&ltogb);
1119:   MatSetLocalToGlobalMapping(J,ltog,ltog);
1120:   MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);
1121: 
1122:   /*
1123:     For each node in the grid: we get the neighbors in the local (on processor ordering
1124:     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1125:     PETSc ordering.
1126:   */
1127:   if (!da->prealloc_only) {
1128:     PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);
1129:     PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));
1130:     for (i=xs; i<xs+nx; i++) {
1131:       istart = PetscMax(-s,gxs - i);
1132:       iend   = PetscMin(s,gxs + gnx - i - 1);
1133:       slot   = i - gxs;
1134: 
1135:       cnt  = 0;
1136:       for (l=0; l<nc; l++) {
1137:         for (i1=istart; i1<iend+1; i1++) {
1138:           cols[cnt++] = l + nc*(slot + i1);
1139:         }
1140:         rows[l]      = l + nc*(slot);
1141:       }
1142:       MatSetValuesLocal(J,nc,rows,cnt,cols,values,INSERT_VALUES);
1143:     }
1144:     PetscFree(values);
1145:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
1146:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
1147:   }
1148:   PetscFree2(rows,cols);
1149:   return(0);
1150: }

1154: PetscErrorCode DMGetMatrix_DA_2d_MPIBAIJ(DM da,Mat J)
1155: {
1156:   PetscErrorCode         ierr;
1157:   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1158:   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
1159:   PetscInt               istart,iend,jstart,jend,ii,jj;
1160:   MPI_Comm               comm;
1161:   PetscScalar            *values;
1162:   DMDABoundaryType       bx,by;
1163:   DMDAStencilType        st;
1164:   ISLocalToGlobalMapping ltog,ltogb;

1167:   /*     
1168:      nc - number of components per grid point 
1169:      col - number of colors needed in one direction for single component problem
1170:   */
1171:   DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);
1172:   col = 2*s + 1;

1174:   DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);
1175:   DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);
1176:   PetscObjectGetComm((PetscObject)da,&comm);

1178:   PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);

1180:   DMGetLocalToGlobalMapping(da,&ltog);
1181:   DMGetLocalToGlobalMappingBlock(da,&ltogb);

1183:   /* determine the matrix preallocation information */
1184:   MatPreallocateInitialize(comm,nx*ny,nx*ny,dnz,onz);
1185:   for (i=xs; i<xs+nx; i++) {
1186:     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1187:     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1188:     for (j=ys; j<ys+ny; j++) {
1189:       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1190:       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1191:       slot = i - gxs + gnx*(j - gys);

1193:       /* Find block columns in block row */
1194:       cnt  = 0;
1195:       for (ii=istart; ii<iend+1; ii++) {
1196:         for (jj=jstart; jj<jend+1; jj++) {
1197:           if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
1198:             cols[cnt++]  = slot + ii + gnx*jj;
1199:           }
1200:         }
1201:       }
1202:       MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);
1203:     }
1204:   }
1205:   MatSeqBAIJSetPreallocation(J,nc,0,dnz);
1206:   MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);
1207:   MatPreallocateFinalize(dnz,onz);

1209:   MatSetLocalToGlobalMapping(J,ltog,ltog);
1210:   MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);

1212:   /*
1213:     For each node in the grid: we get the neighbors in the local (on processor ordering
1214:     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1215:     PETSc ordering.
1216:   */
1217:   if (!da->prealloc_only) {
1218:     PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);
1219:     PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));
1220:     for (i=xs; i<xs+nx; i++) {
1221:       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1222:       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1223:       for (j=ys; j<ys+ny; j++) {
1224:         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1225:         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1226:         slot = i - gxs + gnx*(j - gys);
1227:         cnt  = 0;
1228:         for (ii=istart; ii<iend+1; ii++) {
1229:           for (jj=jstart; jj<jend+1; jj++) {
1230:             if (st == DMDA_STENCIL_BOX || !ii || !jj) { /* BOX or on the STAR */
1231:               cols[cnt++]  = slot + ii + gnx*jj;
1232:             }
1233:           }
1234:         }
1235:         MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);
1236:       }
1237:     }
1238:     PetscFree(values);
1239:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
1240:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
1241:   }
1242:   PetscFree(cols);
1243:   return(0);
1244: }

1248: PetscErrorCode DMGetMatrix_DA_3d_MPIBAIJ(DM da,Mat J)
1249: {
1250:   PetscErrorCode         ierr;
1251:   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1252:   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
1253:   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
1254:   MPI_Comm               comm;
1255:   PetscScalar            *values;
1256:   DMDABoundaryType       bx,by,bz;
1257:   DMDAStencilType        st;
1258:   ISLocalToGlobalMapping ltog,ltogb;

1261:   /*     
1262:          nc - number of components per grid point 
1263:          col - number of colors needed in one direction for single component problem
1264:   
1265:   */
1266:   DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);
1267:   col    = 2*s + 1;

1269:   DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);
1270:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);
1271:   PetscObjectGetComm((PetscObject)da,&comm);

1273:   PetscMalloc(col*col*col*sizeof(PetscInt),&cols);

1275:   DMGetLocalToGlobalMapping(da,&ltog);
1276:   DMGetLocalToGlobalMappingBlock(da,&ltogb);

1278:   /* determine the matrix preallocation information */
1279:   MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);
1280:   for (i=xs; i<xs+nx; i++) {
1281:     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1282:     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1283:     for (j=ys; j<ys+ny; j++) {
1284:       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1285:       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1286:       for (k=zs; k<zs+nz; k++) {
1287:         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1288:         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));

1290:         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);

1292:         /* Find block columns in block row */
1293:         cnt  = 0;
1294:         for (ii=istart; ii<iend+1; ii++) {
1295:           for (jj=jstart; jj<jend+1; jj++) {
1296:             for (kk=kstart; kk<kend+1; kk++) {
1297:               if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1298:                 cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
1299:               }
1300:             }
1301:           }
1302:         }
1303:         MatPreallocateSetLocal(ltogb,1,&slot,ltogb,cnt,cols,dnz,onz);
1304:       }
1305:     }
1306:   }
1307:   MatSeqBAIJSetPreallocation(J,nc,0,dnz);
1308:   MatMPIBAIJSetPreallocation(J,nc,0,dnz,0,onz);
1309:   MatPreallocateFinalize(dnz,onz);

1311:   MatSetLocalToGlobalMapping(J,ltog,ltog);
1312:   MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);

1314:   /*
1315:     For each node in the grid: we get the neighbors in the local (on processor ordering
1316:     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1317:     PETSc ordering.
1318:   */
1319:   if (!da->prealloc_only) {
1320:     PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);
1321:     PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));
1322:     for (i=xs; i<xs+nx; i++) {
1323:       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1324:       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1325:       for (j=ys; j<ys+ny; j++) {
1326:         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1327:         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1328:         for (k=zs; k<zs+nz; k++) {
1329:           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1330:           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
1331: 
1332:           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1333: 
1334:           cnt  = 0;
1335:           for (ii=istart; ii<iend+1; ii++) {
1336:             for (jj=jstart; jj<jend+1; jj++) {
1337:               for (kk=kstart; kk<kend+1; kk++) {
1338:                 if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1339:                   cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
1340:                 }
1341:               }
1342:             }
1343:           }
1344:           MatSetValuesBlockedLocal(J,1,&slot,cnt,cols,values,INSERT_VALUES);
1345:         }
1346:       }
1347:     }
1348:     PetscFree(values);
1349:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
1350:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
1351:   }
1352:   PetscFree(cols);
1353:   return(0);
1354: }

1358: /*
1359:   This helper is for of SBAIJ preallocation, to discard the lower-triangular values which are difficult to
1360:   identify in the local ordering with periodic domain.
1361: */
1362: static PetscErrorCode L2GFilterUpperTriangular(ISLocalToGlobalMapping ltog,PetscInt *row,PetscInt *cnt,PetscInt col[])
1363: {
1365:   PetscInt       i,n;

1368:   ISLocalToGlobalMappingApply(ltog,1,row,row);
1369:   ISLocalToGlobalMappingApply(ltog,*cnt,col,col);
1370:   for (i=0,n=0; i<*cnt; i++) {
1371:     if (col[i] >= *row) col[n++] = col[i];
1372:   }
1373:   *cnt = n;
1374:   return(0);
1375: }

1379: PetscErrorCode DMGetMatrix_DA_2d_MPISBAIJ(DM da,Mat J)
1380: {
1381:   PetscErrorCode         ierr;
1382:   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1383:   PetscInt               m,n,dim,s,*cols,nc,col,cnt,*dnz,*onz;
1384:   PetscInt               istart,iend,jstart,jend,ii,jj;
1385:   MPI_Comm               comm;
1386:   PetscScalar            *values;
1387:   DMDABoundaryType       bx,by;
1388:   DMDAStencilType        st;
1389:   ISLocalToGlobalMapping ltog,ltogb;

1392:   /*     
1393:      nc - number of components per grid point 
1394:      col - number of colors needed in one direction for single component problem
1395:   */
1396:   DMDAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&bx,&by,0,&st);
1397:   col = 2*s + 1;

1399:   DMDAGetCorners(da,&xs,&ys,0,&nx,&ny,0);
1400:   DMDAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);
1401:   PetscObjectGetComm((PetscObject)da,&comm);

1403:   PetscMalloc(col*col*nc*nc*sizeof(PetscInt),&cols);

1405:   DMGetLocalToGlobalMapping(da,&ltog);
1406:   DMGetLocalToGlobalMappingBlock(da,&ltogb);

1408:   /* determine the matrix preallocation information */
1409:   MatPreallocateSymmetricInitialize(comm,nx*ny,nx*ny,dnz,onz);
1410:   for (i=xs; i<xs+nx; i++) {
1411:     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1412:     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1413:     for (j=ys; j<ys+ny; j++) {
1414:       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1415:       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1416:       slot = i - gxs + gnx*(j - gys);

1418:       /* Find block columns in block row */
1419:       cnt  = 0;
1420:       for (ii=istart; ii<iend+1; ii++) {
1421:         for (jj=jstart; jj<jend+1; jj++) {
1422:           if (st == DMDA_STENCIL_BOX || !ii || !jj) {
1423:             cols[cnt++]  = slot + ii + gnx*jj;
1424:           }
1425:         }
1426:       }
1427:       L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);
1428:       MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);
1429:     }
1430:   }
1431:   MatSeqSBAIJSetPreallocation(J,nc,0,dnz);
1432:   MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);
1433:   MatPreallocateFinalize(dnz,onz);

1435:   MatSetLocalToGlobalMapping(J,ltog,ltog);
1436:   MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);

1438:   /*
1439:     For each node in the grid: we get the neighbors in the local (on processor ordering
1440:     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1441:     PETSc ordering.
1442:   */
1443:   if (!da->prealloc_only) {
1444:     PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);
1445:     PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));
1446:     for (i=xs; i<xs+nx; i++) {
1447:       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1448:       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1449:       for (j=ys; j<ys+ny; j++) {
1450:         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1451:         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1452:         slot = i - gxs + gnx*(j - gys);

1454:         /* Find block columns in block row */
1455:         cnt  = 0;
1456:         for (ii=istart; ii<iend+1; ii++) {
1457:           for (jj=jstart; jj<jend+1; jj++) {
1458:             if (st == DMDA_STENCIL_BOX || !ii || !jj) {
1459:               cols[cnt++]  = slot + ii + gnx*jj;
1460:             }
1461:           }
1462:         }
1463:         L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);
1464:         MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);
1465:       }
1466:     }
1467:     PetscFree(values);
1468:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
1469:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
1470:   }
1471:   PetscFree(cols);
1472:   return(0);
1473: }

1477: PetscErrorCode DMGetMatrix_DA_3d_MPISBAIJ(DM da,Mat J)
1478: {
1479:   PetscErrorCode         ierr;
1480:   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1481:   PetscInt               m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
1482:   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
1483:   MPI_Comm               comm;
1484:   PetscScalar            *values;
1485:   DMDABoundaryType       bx,by,bz;
1486:   DMDAStencilType        st;
1487:   ISLocalToGlobalMapping ltog,ltogb;

1490:   /*     
1491:      nc - number of components per grid point 
1492:      col - number of colors needed in one direction for single component problem 
1493:   */
1494:   DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);
1495:   col = 2*s + 1;

1497:   DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);
1498:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);
1499:   PetscObjectGetComm((PetscObject)da,&comm);

1501:   /* create the matrix */
1502:   PetscMalloc(col*col*col*sizeof(PetscInt),&cols);

1504:   DMGetLocalToGlobalMapping(da,&ltog);
1505:   DMGetLocalToGlobalMappingBlock(da,&ltogb);

1507:   /* determine the matrix preallocation information */
1508:   MatPreallocateSymmetricInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);
1509:   for (i=xs; i<xs+nx; i++) {
1510:     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1511:     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1512:     for (j=ys; j<ys+ny; j++) {
1513:       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1514:       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1515:       for (k=zs; k<zs+nz; k++) {
1516:         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1517:         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));

1519:         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);

1521:         /* Find block columns in block row */
1522:         cnt  = 0;
1523:         for (ii=istart; ii<iend+1; ii++) {
1524:           for (jj=jstart; jj<jend+1; jj++) {
1525:             for (kk=kstart; kk<kend+1; kk++) {
1526:               if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
1527:                 cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
1528:               }
1529:             }
1530:           }
1531:         }
1532:         L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);
1533:         MatPreallocateSymmetricSet(slot,cnt,cols,dnz,onz);
1534:       }
1535:     }
1536:   }
1537:   MatSeqSBAIJSetPreallocation(J,nc,0,dnz);
1538:   MatMPISBAIJSetPreallocation(J,nc,0,dnz,0,onz);
1539:   MatPreallocateFinalize(dnz,onz);

1541:   MatSetLocalToGlobalMapping(J,ltog,ltog);
1542:   MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);

1544:   /*
1545:     For each node in the grid: we get the neighbors in the local (on processor ordering
1546:     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1547:     PETSc ordering.
1548:   */
1549:   if (!da->prealloc_only) {
1550:     PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);
1551:     PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));
1552:     for (i=xs; i<xs+nx; i++) {
1553:       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1554:       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1555:       for (j=ys; j<ys+ny; j++) {
1556:         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1557:         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1558:         for (k=zs; k<zs+nz; k++) {
1559:           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1560:           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
1561: 
1562:           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1563: 
1564:           cnt  = 0;
1565:           for (ii=istart; ii<iend+1; ii++) {
1566:             for (jj=jstart; jj<jend+1; jj++) {
1567:               for (kk=kstart; kk<kend+1; kk++) {
1568:                 if ((st == DMDA_STENCIL_BOX) || (!ii && !jj) || (!jj && !kk) || (!ii && !kk)) {
1569:                   cols[cnt++]  = slot + ii + gnx*jj + gnx*gny*kk;
1570:                 }
1571:               }
1572:             }
1573:           }
1574:           L2GFilterUpperTriangular(ltogb,&slot,&cnt,cols);
1575:           MatSetValuesBlocked(J,1,&slot,cnt,cols,values,INSERT_VALUES);
1576:         }
1577:       }
1578:     }
1579:     PetscFree(values);
1580:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
1581:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
1582:   }
1583:   PetscFree(cols);
1584:   return(0);
1585: }

1587: /* ---------------------------------------------------------------------------------*/

1591: PetscErrorCode DMGetMatrix_DA_3d_MPIAIJ_Fill(DM da,Mat J)
1592: {
1593:   PetscErrorCode         ierr;
1594:   PetscInt               xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
1595:   PetscInt               m,n,dim,s,*cols,k,nc,row,col,cnt,l,p,*dnz,*onz;
1596:   PetscInt               istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
1597:   DM_DA                  *dd = (DM_DA*)da->data;
1598:   PetscInt               ifill_col,*dfill = dd->dfill,*ofill = dd->ofill;
1599:   MPI_Comm               comm;
1600:   PetscScalar            *values;
1601:   DMDABoundaryType       bx,by,bz;
1602:   ISLocalToGlobalMapping ltog,ltogb;
1603:   DMDAStencilType        st;

1606:   /*     
1607:          nc - number of components per grid point 
1608:          col - number of colors needed in one direction for single component problem
1609:   
1610:   */
1611:   DMDAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&bx,&by,&bz,&st);
1612:   col    = 2*s + 1;
1613:   if (bx == DMDA_BOUNDARY_PERIODIC && (m % col)){
1614:     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisible\n\
1615:                  by 2*stencil_width + 1\n");
1616:   }
1617:   if (by == DMDA_BOUNDARY_PERIODIC && (n % col)){
1618:     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisible\n\
1619:                  by 2*stencil_width + 1\n");
1620:   }
1621:   if (bz == DMDA_BOUNDARY_PERIODIC && (p % col)){
1622:     SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisible\n\
1623:                  by 2*stencil_width + 1\n");
1624:   }

1626:   DMDAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);
1627:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);
1628:   PetscObjectGetComm((PetscObject)da,&comm);

1630:   PetscMalloc(col*col*col*nc*sizeof(PetscInt),&cols);
1631:   DMGetLocalToGlobalMapping(da,&ltog);
1632:   DMGetLocalToGlobalMappingBlock(da,&ltogb);

1634:   /* determine the matrix preallocation information */
1635:   MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);


1638:   for (i=xs; i<xs+nx; i++) {
1639:     istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1640:     iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1641:     for (j=ys; j<ys+ny; j++) {
1642:       jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1643:       jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1644:       for (k=zs; k<zs+nz; k++) {
1645:         kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1646:         kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
1647: 
1648:         slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1649: 
1650:         for (l=0; l<nc; l++) {
1651:           cnt  = 0;
1652:           for (ii=istart; ii<iend+1; ii++) {
1653:             for (jj=jstart; jj<jend+1; jj++) {
1654:               for (kk=kstart; kk<kend+1; kk++) {
1655:                 if (ii || jj || kk) {
1656:                   if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1657:                     for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
1658:                       cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1659:                   }
1660:                 } else {
1661:                   if (dfill) {
1662:                     for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
1663:                       cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1664:                   } else {
1665:                     for (ifill_col=0; ifill_col<nc; ifill_col++)
1666:                       cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1667:                   }
1668:                 }
1669:               }
1670:             }
1671:           }
1672:           row  = l + nc*(slot);
1673:           MatPreallocateSetLocal(ltog,1,&row,ltog,cnt,cols,dnz,onz);
1674:         }
1675:       }
1676:     }
1677:   }
1678:   MatSeqAIJSetPreallocation(J,0,dnz);
1679:   MatMPIAIJSetPreallocation(J,0,dnz,0,onz);
1680:   MatPreallocateFinalize(dnz,onz);
1681:   MatSetLocalToGlobalMapping(J,ltog,ltog);
1682:   MatSetLocalToGlobalMappingBlock(J,ltogb,ltogb);

1684:   /*
1685:     For each node in the grid: we get the neighbors in the local (on processor ordering
1686:     that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
1687:     PETSc ordering.
1688:   */
1689:   if (!da->prealloc_only) {
1690:     PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);
1691:     PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));
1692:     for (i=xs; i<xs+nx; i++) {
1693:       istart = (bx == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-i));
1694:       iend   = (bx == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,m-i-1));
1695:       for (j=ys; j<ys+ny; j++) {
1696:         jstart = (by == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-j));
1697:         jend   = (by == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,n-j-1));
1698:         for (k=zs; k<zs+nz; k++) {
1699:           kstart = (bz == DMDA_BOUNDARY_PERIODIC) ? -s : (PetscMax(-s,-k));
1700:           kend   = (bz == DMDA_BOUNDARY_PERIODIC) ?  s : (PetscMin(s,p-k-1));
1701: 
1702:           slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
1703: 
1704:           for (l=0; l<nc; l++) {
1705:             cnt  = 0;
1706:             for (ii=istart; ii<iend+1; ii++) {
1707:               for (jj=jstart; jj<jend+1; jj++) {
1708:                 for (kk=kstart; kk<kend+1; kk++) {
1709:                   if (ii || jj || kk) {
1710:                     if ((st == DMDA_STENCIL_BOX) || ((!ii && !jj) || (!jj && !kk) || (!ii && !kk))) {/* entries on star*/
1711:                       for (ifill_col=ofill[l]; ifill_col<ofill[l+1]; ifill_col++)
1712:                         cols[cnt++]  = ofill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1713:                     }
1714:                   } else {
1715:                     if (dfill) {
1716:                       for (ifill_col=dfill[l]; ifill_col<dfill[l+1]; ifill_col++)
1717:                         cols[cnt++]  = dfill[ifill_col] + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1718:                     } else {
1719:                       for (ifill_col=0; ifill_col<nc; ifill_col++)
1720:                         cols[cnt++]  = ifill_col + nc*(slot + ii + gnx*jj + gnx*gny*kk);
1721:                     }
1722:                   }
1723:                 }
1724:               }
1725:             }
1726:             row  = l + nc*(slot);
1727:             MatSetValuesLocal(J,1,&row,cnt,cols,values,INSERT_VALUES);
1728:           }
1729:         }
1730:       }
1731:     }
1732:     PetscFree(values);
1733:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
1734:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
1735:   }
1736:   PetscFree(cols);
1737:   return(0);
1738: }