Actual source code: fdda.c
1: /*$Id: fdda.c,v 1.65 2001/04/10 19:37:31 bsmith Exp $*/
2:
3: #include "petscda.h" /*I "petscda.h" I*/
4: #include "petscmat.h" /*I "petscmat.h" I*/
5: #include "src/dm/da/daimpl.h"
7: EXTERN int DAGetColoring1d_MPIAIJ(DA,ISColoringType,ISColoring *,Mat *);
8: EXTERN int DAGetColoring2d_MPIAIJ(DA,ISColoringType,ISColoring *,Mat *);
9: EXTERN int DAGetColoring3d_MPIAIJ(DA,ISColoringType,ISColoring *,Mat *);
10: EXTERN int DAGetColoring3d_MPIBAIJ(DA,ISColoringType,ISColoring *,Mat *);
12: /*@C
13: DAGetColoring - Gets the coloring and nonzero structure required for computing the Jacobian via
14: finite differences on a function defined using a stencil on the DA.
16: Collective on DA
18: Input Parameter:
19: + da - the distributed array
20: - mtype - either MATMPIAIJ or MATMPIBAIJ
22: Output Parameters:
23: + coloring - matrix coloring for use in computing Jacobians (or PETSC_NULL if not needed)
24: - J - matrix with the correct nonzero structure (or PETSC_NULL if not needed)
25: (obviously without the correct Jacobian values)
27: Level: advanced
29: Notes: These compute the graph coloring of the graph of A^{T}A. The coloring used
30: for efficient (parallel or thread based) triangular solves etc is NOT yet
31: available.
33: This does not yet handle BAIJ matrices, because
34: 1) we need a way for the user to indicate what matrix type they want
35: 2) after the matrix is created, for BAIJ matrices we need to set nc to 1 and
36: use MatSetValuesBlockedLocal() instead of MatSetValuesLocal()
38: .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), DAGetColoringMPIBAIJ()
40: @*/
41: int DAGetColoring(DA da,ISColoringType ctype,MatType mtype,ISColoring *coloring,Mat *J)
42: {
43: int ierr,dim;
44: PetscTruth aij;
47: /*
48: m
49: ------------------------------------------------------
50: | |
51: | |
52: | ---------------------- |
53: | | | |
54: n | yn | | |
55: | | | |
56: | .--------------------- |
57: | (xs,ys) xn |
58: | . |
59: | (gxs,gys) |
60: | |
61: -----------------------------------------------------
62: */
64: /*
65: nc - number of components per grid point
66: col - number of colors needed in one direction for single component problem
67:
68: */
69: DAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0);
71: /*
72: We do not provide a getcoloring function in the DA operations because
73: the basic DA does not know about matrices. We think of DA as being more
74: more low-level then matrices.
75: */
76: PetscStrcmp(MATMPIAIJ,mtype,&aij);
77: if (aij) {
78: if (dim == 1) {
79: DAGetColoring1d_MPIAIJ(da,ctype,coloring,J);
80: } else if (dim == 2) {
81: DAGetColoring2d_MPIAIJ(da,ctype,coloring,J);
82: } else if (dim == 3) {
83: DAGetColoring3d_MPIAIJ(da,ctype,coloring,J);
84: }
85: } else {
86: if (dim == 3) {
87: DAGetColoring3d_MPIBAIJ(da,ctype,coloring,J);
88: } else {
89: SETERRQ1(1,"Not done for %d dimension, send use mail petsc-maint@mcs.anl.gov for code",dim);
90: }
91: }
92: return(0);
93: }
95: /* ---------------------------------------------------------------------------------*/
97: int DAGetColoring2d_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring,Mat *J)
98: {
99: int ierr,xs,ys,nx,ny,*colors,i,j,ii,slot,gxs,gys,gnx,gny;
100: int m,n,dim,w,s,*cols,k,nc,*rows,col,cnt,l,p;
101: int lstart,lend,pstart,pend,*dnz,*onz,size;
102: MPI_Comm comm;
103: Scalar *values;
104: DAPeriodicType wrap;
105: ISLocalToGlobalMapping ltog;
106: DAStencilType st;
109: /*
110: nc - number of components per grid point
111: col - number of colors needed in one direction for single component problem
112:
113: */
114: DAGetInfo(da,&dim,&m,&n,0,0,0,0,&w,&s,&wrap,&st);
115: nc = w;
116: col = 2*s + 1;
117: if (DAXPeriodic(wrap) && (m % col)){
118: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisiblen
119: by 2*stencil_width + 1n");
120: }
121: if (DAYPeriodic(wrap) && (n % col)){
122: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisiblen
123: by 2*stencil_width + 1n");
124: }
125: DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);
126: DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);
127: PetscObjectGetComm((PetscObject)da,&comm);
128: MPI_Comm_size(comm,&size);
130: /* create the coloring */
131: if (coloring) {
132: if (ctype == IS_COLORING_GLOBAL) {
133: PetscMalloc(nc*nx*ny*sizeof(int),&colors);
134: ii = 0;
135: for (j=ys; j<ys+ny; j++) {
136: for (i=xs; i<xs+nx; i++) {
137: for (k=0; k<nc; k++) {
138: colors[ii++] = k + nc*((i % col) + col*(j % col));
139: }
140: }
141: }
142: ISColoringCreate(comm,nc*nx*ny,colors,coloring);
143: } else if (ctype == IS_COLORING_LOCAL) {
144: PetscMalloc(nc*gnx*gny*sizeof(int),&colors);
145: ii = 0;
146: for (j=gys; j<gys+gny; j++) {
147: for (i=gxs; i<gxs+gnx; i++) {
148: for (k=0; k<nc; k++) {
149: colors[ii++] = k + nc*((i % col) + col*(j % col));
150: }
151: }
152: }
153: ISColoringCreate(comm,nc*gnx*gny,colors,coloring);
154: } else SETERRQ1(1,"Unknown ISColoringType %d",ctype);
155: }
157: /* Create the matrix */
158: if (J) {
159: int bs = nc,dims[2],starts[2];
160: /* create empty Jacobian matrix */
161: ierr = MatCreate(comm,nc*nx*ny,nc*nx*ny,PETSC_DECIDE,PETSC_DECIDE,J);
163: PetscMalloc(col*col*nc*nc*sizeof(Scalar),&values);
164: PetscMemzero(values,col*col*nc*nc*sizeof(Scalar));
165: PetscMalloc(nc*sizeof(int),&rows);
166: PetscMalloc(col*col*nc*nc*sizeof(int),&cols);
167: DAGetISLocalToGlobalMapping(da,<og);
169: /* determine the matrix preallocation information */
170: MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);
171: for (i=xs; i<xs+nx; i++) {
173: pstart = PetscMax(-s,-i);
174: pend = PetscMin(s,m-i-1);
176: for (j=ys; j<ys+ny; j++) {
177: slot = i - gxs + gnx*(j - gys);
179: lstart = PetscMax(-s,-j);
180: lend = PetscMin(s,n-j-1);
182: cnt = 0;
183: for (k=0; k<nc; k++) {
184: for (l=lstart; l<lend+1; l++) {
185: for (p=pstart; p<pend+1; p++) {
186: if ((st == DA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */
187: cols[cnt++] = k + nc*(slot + gnx*l + p);
188: }
189: }
190: }
191: rows[k] = k + nc*(slot);
192: }
193: MatPreallocateSetLocal(ltog,nc,rows,cnt,cols,dnz,onz);
194: }
195: }
196: /* set matrix type and preallocation information */
197: if (size > 1) {
198: MatSetType(*J,MATMPIAIJ);
199: } else {
200: MatSetType(*J,MATSEQAIJ);
201: }
202: MatSeqAIJSetPreallocation(*J,0,dnz);
203: MatSeqBAIJSetPreallocation(*J,bs,0,dnz);
204: MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);
205: MatMPIBAIJSetPreallocation(*J,bs,0,dnz,0,onz);
206: MatPreallocateFinalize(dnz,onz);
207: MatSetLocalToGlobalMapping(*J,ltog);
208: DAGetGhostCorners(da,&starts[0],&starts[1],PETSC_IGNORE,&dims[0],&dims[1],PETSC_IGNORE);
209: MatSetStencil(*J,2,dims,starts,nc);
211: /*
212: For each node in the grid: we get the neighbors in the local (on processor ordering
213: that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
214: PETSc ordering.
215: */
216: for (i=xs; i<xs+nx; i++) {
218: pstart = PetscMax(-s,-i);
219: pend = PetscMin(s,m-i-1);
221: for (j=ys; j<ys+ny; j++) {
222: slot = i - gxs + gnx*(j - gys);
224: lstart = PetscMax(-s,-j);
225: lend = PetscMin(s,n-j-1);
227: cnt = 0;
228: for (k=0; k<nc; k++) {
229: for (l=lstart; l<lend+1; l++) {
230: for (p=pstart; p<pend+1; p++) {
231: if ((st == DA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */
232: cols[cnt++] = k + nc*(slot + gnx*l + p);
233: }
234: }
235: }
236: rows[k] = k + nc*(slot);
237: }
238: MatSetValuesLocal(*J,nc,rows,cnt,cols,values,INSERT_VALUES);
239: }
240: }
241: PetscFree(values);
242: PetscFree(rows);
243: PetscFree(cols);
244: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
245: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
246: }
247: return(0);
248: }
250: /* ---------------------------------------------------------------------------------*/
252: int DAGetColoring3d_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring,Mat *J)
253: {
254: int ierr,xs,ys,nx,ny,*colors,i,j,slot,gxs,gys,gnx,gny;
255: int m,n,dim,s,*cols,k,nc,*rows,col,cnt,l,p,*dnz,*onz;
256: int istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk,size;
257: MPI_Comm comm;
258: Scalar *values;
259: DAPeriodicType wrap;
260: ISLocalToGlobalMapping ltog;
261: DAStencilType st;
264: /*
265: nc - number of components per grid point
266: col - number of colors needed in one direction for single component problem
267:
268: */
269: DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);
270: col = 2*s + 1;
271: if (DAXPeriodic(wrap) && (m % col)){
272: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisiblen
273: by 2*stencil_width + 1n");
274: }
275: if (DAYPeriodic(wrap) && (n % col)){
276: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisiblen
277: by 2*stencil_width + 1n");
278: }
279: if (DAZPeriodic(wrap) && (p % col)){
280: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisiblen
281: by 2*stencil_width + 1n");
282: }
284: DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);
285: DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);
286: PetscObjectGetComm((PetscObject)da,&comm);
287: MPI_Comm_size(comm,&size);
289: /* create the coloring */
290: if (coloring) {
291: if (ctype == IS_COLORING_GLOBAL) {
292: PetscMalloc(nc*nx*ny*nz*sizeof(int),&colors);
293: ii = 0;
294: for (k=zs; k<zs+nz; k++) {
295: for (j=ys; j<ys+ny; j++) {
296: for (i=xs; i<xs+nx; i++) {
297: for (l=0; l<nc; l++) {
298: colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
299: }
300: }
301: }
302: }
303: ISColoringCreate(comm,nc*nx*ny*nz,colors,coloring);
304: } else if (ctype == IS_COLORING_LOCAL) {
305: PetscMalloc(nc*gnx*gny*gnz*sizeof(int),&colors);
306: ii = 0;
307: for (k=gzs; k<gzs+gnz; k++) {
308: for (j=gys; j<gys+gny; j++) {
309: for (i=gxs; i<gxs+gnx; i++) {
310: for (l=0; l<nc; l++) {
311: colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
312: }
313: }
314: }
315: }
316: ISColoringCreate(comm,nc*gnx*gny*gnz,colors,coloring);
317: } else SETERRQ1(1,"Unknown ISColoringType %d",ctype);
318: }
320: /* create the matrix */
321: if (J) {
322: int bs = nc,dims[3],starts[3];
323: /* create empty Jacobian matrix */
324: MatCreate(comm,nc*nx*ny*nz,nc*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE,J);
325: PetscMalloc(col*col*col*nc*nc*nc*sizeof(Scalar),&values);
326: PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(Scalar));
327: PetscMalloc(nc*sizeof(int),&rows);
328: PetscMalloc(col*col*col*nc*sizeof(int),&cols);
329: DAGetISLocalToGlobalMapping(da,<og);
331: /* determine the matrix preallocation information */
332: MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);
333: for (i=xs; i<xs+nx; i++) {
334: istart = PetscMax(-s,-i);
335: iend = PetscMin(s,m-i-1);
336: for (j=ys; j<ys+ny; j++) {
337: jstart = PetscMax(-s,-j);
338: jend = PetscMin(s,n-j-1);
339: for (k=zs; k<zs+nz; k++) {
340: kstart = PetscMax(-s,-k);
341: kend = PetscMin(s,p-k-1);
343: slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
345: cnt = 0;
346: for (l=0; l<nc; l++) {
347: for (ii=istart; ii<iend+1; ii++) {
348: for (jj=jstart; jj<jend+1; jj++) {
349: for (kk=kstart; kk<kend+1; kk++) {
350: if ((st == DA_STENCIL_BOX) || (!ii || !jj || !kk)) { /* entries on star */
351: cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
352: }
353: }
354: }
355: }
356: rows[l] = l + nc*(slot);
357: }
358: MatPreallocateSetLocal(ltog,nc,rows,cnt,cols,dnz,onz);
359: }
360: }
361: }
362: /* set matrix type and preallocation */
363: if (size > 1) {
364: MatSetType(*J,MATMPIAIJ);
365: } else {
366: MatSetType(*J,MATSEQAIJ);
367: }
368: MatSeqAIJSetPreallocation(*J,0,dnz);
369: MatSeqBAIJSetPreallocation(*J,bs,0,dnz);
370: MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);
371: MatMPIBAIJSetPreallocation(*J,bs,0,dnz,0,onz);
372: MatPreallocateFinalize(dnz,onz);
373: MatSetLocalToGlobalMapping(*J,ltog);
374: DAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);
375: MatSetStencil(*J,3,dims,starts,nc);
377: /*
378: For each node in the grid: we get the neighbors in the local (on processor ordering
379: that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
380: PETSc ordering.
381: */
382: for (i=xs; i<xs+nx; i++) {
383: istart = PetscMax(-s,-i);
384: iend = PetscMin(s,m-i-1);
385: for (j=ys; j<ys+ny; j++) {
386: jstart = PetscMax(-s,-j);
387: jend = PetscMin(s,n-j-1);
388: for (k=zs; k<zs+nz; k++) {
389: kstart = PetscMax(-s,-k);
390: kend = PetscMin(s,p-k-1);
392: slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
394: cnt = 0;
395: for (l=0; l<nc; l++) {
396: for (ii=istart; ii<iend+1; ii++) {
397: for (jj=jstart; jj<jend+1; jj++) {
398: for (kk=kstart; kk<kend+1; kk++) {
399: if ((st == DA_STENCIL_BOX) || (!ii || !jj || !kk)) { /* entries on star */
400: cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
401: }
402: }
403: }
404: }
405: rows[l] = l + nc*(slot);
406: }
407: MatSetValuesLocal(*J,nc,rows,cnt,cols,values,INSERT_VALUES);
408: }
409: }
410: }
411: PetscFree(values);
412: PetscFree(rows);
413: PetscFree(cols);
414: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
415: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
416: }
417: return(0);
418: }
420: /* ---------------------------------------------------------------------------------*/
422: int DAGetColoring1d_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring,Mat *J)
423: {
424: int ierr,xs,nx,*colors,i,i1,slot,gxs,gnx;
425: int m,dim,s,*cols,nc,*rows,col,cnt,l;
426: int istart,iend,size;
427: MPI_Comm comm;
428: Scalar *values;
429: DAPeriodicType wrap;
430: ISLocalToGlobalMapping ltog;
433: /*
434: nc - number of components per grid point
435: col - number of colors needed in one direction for single component problem
436:
437: */
438: DAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&wrap,0);
439: col = 2*s + 1;
441: if (DAXPeriodic(wrap) && (m % col)) {
442: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points is divisiblen
443: by 2*stencil_width + 1n");
444: }
447: DAGetCorners(da,&xs,0,0,&nx,0,0);
448: DAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);
449: PetscObjectGetComm((PetscObject)da,&comm);
450: MPI_Comm_size(comm,&size);
452: /* create the coloring */
453: if (coloring) {
454: if (ctype == IS_COLORING_GLOBAL) {
455: PetscMalloc(nc*nx*sizeof(int),&colors);
456: i1 = 0;
457: for (i=xs; i<xs+nx; i++) {
458: for (l=0; l<nc; l++) {
459: colors[i1++] = l + nc*(i % col);
460: }
461: }
462: ISColoringCreate(comm,nc*nx,colors,coloring);
463: } else if (ctype == IS_COLORING_LOCAL) {
464: PetscMalloc(nc*gnx*sizeof(int),&colors);
465: i1 = 0;
466: for (i=gxs; i<gxs+gnx; i++) {
467: for (l=0; l<nc; l++) {
468: colors[i1++] = l + nc*(i % col);
469: }
470: }
471: ISColoringCreate(comm,nc*gnx,colors,coloring);
472: } else SETERRQ1(1,"Unknown ISColoringType %d",ctype);
473: }
475: /* create empty Jacobian matrix */
476: if (J) {
477: int bs = nc,dims[1],starts[1];
478:
479: ierr = MatCreate(comm,nc*nx,nc*nx,PETSC_DECIDE,PETSC_DECIDE,J);
480: if (size > 1) {
481: MatSetType(*J,MATMPIAIJ);
482: } else {
483: MatSetType(*J,MATSEQAIJ);
484: }
485: MatSeqAIJSetPreallocation(*J,col*nc,0);
486: MatSeqBAIJSetPreallocation(*J,bs,col,0);
487: MatMPIAIJSetPreallocation(*J,col*nc,0,0,0);
488: MatMPIBAIJSetPreallocation(*J,bs,col,0,0,0);
489: DAGetGhostCorners(da,&starts[0],PETSC_IGNORE,PETSC_IGNORE,&dims[0],PETSC_IGNORE,PETSC_IGNORE);
490: MatSetStencil(*J,1,dims,starts,nc);
492: PetscMalloc(col*nc*nc*sizeof(Scalar),&values);
493: PetscMemzero(values,col*nc*nc*sizeof(Scalar));
494: PetscMalloc(nc*sizeof(int),&rows);
495: PetscMalloc(col*nc*sizeof(int),&cols);
496:
497: DAGetISLocalToGlobalMapping(da,<og);
498: MatSetLocalToGlobalMapping(*J,ltog);
500: /*
501: For each node in the grid: we get the neighbors in the local (on processor ordering
502: that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
503: PETSc ordering.
504: */
505: for (i=xs; i<xs+nx; i++) {
506: istart = PetscMax(-s,gxs - i);
507: iend = PetscMin(s,gxs + gnx - i - 1);
508: slot = i - gxs;
510: cnt = 0;
511: for (l=0; l<nc; l++) {
512: for (i1=istart; i1<iend+1; i1++) {
513: cols[cnt++] = l + nc*(slot + i1);
514: }
515: rows[l] = l + nc*(slot);
516: }
517: MatSetValuesLocal(*J,nc,rows,cnt,cols,values,INSERT_VALUES);
518: }
519: PetscFree(values);
520: PetscFree(rows);
521: PetscFree(cols);
522: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
523: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
524: }
525: return(0);
526: }
528: int DAGetColoring3d_MPIBAIJ(DA da,ISColoringType ctype,ISColoring *coloring,Mat *J)
529: {
530: int ierr,xs,ys,nx,ny,*colors,i,j,slot,gxs,gys,gnx,gny;
531: int m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
532: int istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
533: MPI_Comm comm;
534: Scalar *values;
535: DAPeriodicType wrap;
536: DAStencilType st;
537: ISLocalToGlobalMapping ltog;
540: /*
541: nc - number of components per grid point
542: col - number of colors needed in one direction for single component problem
543:
544: */
545: DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);
546: if (wrap != DA_NONPERIODIC) SETERRQ(PETSC_ERR_SUP,"Currently no support for periodic");
547: col = 2*s + 1;
549: DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);
550: DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);
551: PetscObjectGetComm((PetscObject)da,&comm);
553: /* create the coloring */
554: if (coloring) {
555: if (ctype == IS_COLORING_GLOBAL) {
556: PetscMalloc(nx*ny*nz*sizeof(int),&colors);
557: ii = 0;
558: for (k=zs; k<zs+nz; k++) {
559: for (j=ys; j<ys+ny; j++) {
560: for (i=xs; i<xs+nx; i++) {
561: colors[ii++] = (i % col) + col*(j % col) + col*col*(k % col);
562: }
563: }
564: }
565: ISColoringCreate(comm,nx*ny*nz,colors,coloring);
566: } else if (ctype == IS_COLORING_LOCAL) {
567: PetscMalloc(gnx*gny*gnz*sizeof(int),&colors);
568: ii = 0;
569: for (k=gzs; k<gzs+gnz; k++) {
570: for (j=gys; j<gys+gny; j++) {
571: for (i=gxs; i<gxs+gnx; i++) {
572: colors[ii++] = (i % col) + col*(j % col) + col*col*(k % col);
573: }
574: }
575: }
576: ISColoringCreate(comm,gnx*gny*gnz,colors,coloring);
577: } else SETERRQ1(1,"Unknown ISColoringType %d",ctype);
578: }
580: /* create the matrix */
581: if (J) {
583: ierr = PetscMalloc(col*col*col*nc*nc*sizeof(Scalar),&values);
584: ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(Scalar));
585: ierr = PetscMalloc(col*col*col*sizeof(int),&cols);
587: DAGetISLocalToGlobalMappingBlck(da,<og);
589: /* determine the matrix preallocation information */
590: MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);
591: for (i=xs; i<xs+nx; i++) {
592: istart = PetscMax(-s,-i);
593: iend = PetscMin(s,m-i-1);
594: for (j=ys; j<ys+ny; j++) {
595: jstart = PetscMax(-s,-j);
596: jend = PetscMin(s,n-j-1);
597: for (k=zs; k<zs+nz; k++) {
598: kstart = PetscMax(-s,-k);
599: kend = PetscMin(s,p-k-1);
601: slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
603: /* Find block columns in block row */
604: cnt = 0;
605: if (st == DA_STENCIL_BOX) { /* if using BOX stencil */
606: for (ii=istart; ii<iend+1; ii++) {
607: for (jj=jstart; jj<jend+1; jj++) {
608: for (kk=kstart; kk<kend+1; kk++) {
609: cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
610: }
611: }
612: }
613: } else { /* Star stencil */
614: cnt = 0;
615: for (ii=istart; ii<iend+1; ii++) {
616: if (ii) {
617: /* jj and kk must be zero */
618: /* cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; */
619: cols[cnt++] = slot + ii;
620: } else {
621: for (jj=jstart; jj<jend+1; jj++) {
622: if (jj) {
623: /* ii and kk must be zero */
624: cols[cnt++] = slot + gnx*jj;
625: } else {
626: /* ii and jj must be zero */
627: for (kk=kstart; kk<kend+1; kk++) {
628: cols[cnt++] = slot + gnx*gny*kk;
629: }
630: }
631: }
632: }
633: }
634: }
635: MatPreallocateSetLocal(ltog,1,&slot,cnt,cols,dnz,onz);
636: }
637: }
638: }
640: /* create empty Jacobian matrix */
641: MatCreateMPIBAIJ(comm,nc,nc*nx*ny*nz,nc*nx*ny*nz,PETSC_DECIDE,
642: PETSC_DECIDE,0,dnz,0,onz,J);
644: MatPreallocateFinalize(dnz,onz);
645: MatSetLocalToGlobalMappingBlock(*J,ltog);
647: /*
648: For each node in the grid: we get the neighbors in the local (on processor ordering
649: that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
650: PETSc ordering.
651: */
653: for (i=xs; i<xs+nx; i++) {
654: istart = PetscMax(-s,-i);
655: iend = PetscMin(s,m-i-1);
656: for (j=ys; j<ys+ny; j++) {
657: jstart = PetscMax(-s,-j);
658: jend = PetscMin(s,n-j-1);
659: for (k=zs; k<zs+nz; k++) {
660: kstart = PetscMax(-s,-k);
661: kend = PetscMin(s,p-k-1);
663: slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
665: cnt = 0;
666: if (st == DA_STENCIL_BOX) { /* if using BOX stencil */
667: for (ii=istart; ii<iend+1; ii++) {
668: for (jj=jstart; jj<jend+1; jj++) {
669: for (kk=kstart; kk<kend+1; kk++) {
670: cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
671: }
672: }
673: }
674: } else { /* Star stencil */
675: cnt = 0;
676: for (ii=istart; ii<iend+1; ii++) {
677: if (ii) {
678: /* jj and kk must be zero */
679: /* cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; */
680: cols[cnt++] = slot + ii;
681: } else {
682: for (jj=jstart; jj<jend+1; jj++) {
683: if (jj) {
684: /* ii and kk must be zero */
685: cols[cnt++] = slot + gnx*jj;
686: } else {
687: /* ii and jj must be zero */
688: for (kk=kstart; kk<kend+1; kk++) {
689: cols[cnt++] = slot + gnx*gny*kk;
690: }
691: }
692: }
693: }
694: }
695: }
696: MatSetValuesBlockedLocal(*J,1,&slot,cnt,cols,values,INSERT_VALUES);
697: }
698: }
699: }
700: PetscFree(values);
701: PetscFree(cols);
702: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
703: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
704: }
705: return(0);
706: }