Actual source code: fdda.c
1: /*$Id: fdda.c,v 1.75 2001/08/07 21:31:51 bsmith Exp $*/
2:
3: #include src/dm/da/daimpl.h
4: #include petscmat.h
7: EXTERN int DAGetColoring1d_MPIAIJ(DA,ISColoringType,ISColoring *);
8: EXTERN int DAGetColoring2d_MPIAIJ(DA,ISColoringType,ISColoring *);
9: EXTERN int DAGetColoring2d_5pt_MPIAIJ(DA,ISColoringType,ISColoring *);
10: EXTERN int DAGetColoring3d_MPIAIJ(DA,ISColoringType,ISColoring *);
11: EXTERN int DAGetColoring3d_MPIBAIJ(DA,ISColoringType,ISColoring *);
13: /*@C
14: DAGetColoring - Gets the coloring required for computing the Jacobian via
15: finite differences on a function defined using a stencil on the DA.
17: Collective on DA
19: Input Parameter:
20: + da - the distributed array
21: - ctype - IS_COLORING_LOCAL or IS_COLORING_GHOSTED
23: Output Parameters:
24: . coloring - matrix coloring for use in computing Jacobians (or PETSC_NULL if not needed)
26: Level: advanced
28: Notes: These compute the graph coloring of the graph of A^{T}A. The coloring used
29: for efficient (parallel or thread based) triangular solves etc is NOT yet
30: available.
33: .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), ISColoringType, ISColoring
35: @*/
36: int DAGetColoring(DA da,ISColoringType ctype,ISColoring *coloring)
37: {
38: int ierr,dim;
41: /*
42: m
43: ------------------------------------------------------
44: | |
45: | |
46: | ---------------------- |
47: | | | |
48: n | yn | | |
49: | | | |
50: | .--------------------- |
51: | (xs,ys) xn |
52: | . |
53: | (gxs,gys) |
54: | |
55: -----------------------------------------------------
56: */
58: /*
59: nc - number of components per grid point
60: col - number of colors needed in one direction for single component problem
61:
62: */
63: DAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0);
65: /*
66: We do not provide a getcoloring function in the DA operations because
67: the basic DA does not know about matrices. We think of DA as being more
68: more low-level then matrices.
69: */
70: if (dim == 1) {
71: DAGetColoring1d_MPIAIJ(da,ctype,coloring);
72: } else if (dim == 2) {
73: DAGetColoring2d_MPIAIJ(da,ctype,coloring);
74: } else if (dim == 3) {
75: DAGetColoring3d_MPIAIJ(da,ctype,coloring);
76: } else {
77: SETERRQ1(1,"Not done for %d dimension, send us mail petsc-maint@mcs.anl.gov for code",dim);
78: }
79: return(0);
80: }
82: /* ---------------------------------------------------------------------------------*/
84: int DAGetColoring2d_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring)
85: {
86: int ierr,xs,ys,nx,ny,*colors,i,j,ii,gxs,gys,gnx,gny;
87: int m,n,dim,w,s,k,nc,col,size;
88: MPI_Comm comm;
89: DAPeriodicType wrap;
90: DAStencilType st;
93: /*
94: nc - number of components per grid point
95: col - number of colors needed in one direction for single component problem
96:
97: */
98: DAGetInfo(da,&dim,&m,&n,0,0,0,0,&w,&s,&wrap,&st);
99: nc = w;
100: col = 2*s + 1;
101: if (DAXPeriodic(wrap) && (m % col)){
102: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisiblen
103: by 2*stencil_width + 1n");
104: }
105: if (DAYPeriodic(wrap) && (n % col)){
106: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisiblen
107: by 2*stencil_width + 1n");
108: }
109: DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);
110: DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);
111: PetscObjectGetComm((PetscObject)da,&comm);
112: MPI_Comm_size(comm,&size);
114: /* create the coloring */
115: if (st == DA_STENCIL_STAR && s == 1) {
116: DAGetColoring2d_5pt_MPIAIJ(da,ctype,coloring);
117: } else if (ctype == IS_COLORING_LOCAL) {
118: if (!da->localcoloring) {
119: PetscMalloc(nc*nx*ny*sizeof(int),&colors);
120: ii = 0;
121: for (j=ys; j<ys+ny; j++) {
122: for (i=xs; i<xs+nx; i++) {
123: for (k=0; k<nc; k++) {
124: colors[ii++] = k + nc*((i % col) + col*(j % col));
125: }
126: }
127: }
128: ISColoringCreate(comm,nc*nx*ny,colors,&da->localcoloring);
129: }
130: *coloring = da->localcoloring;
131: } else if (ctype == IS_COLORING_GHOSTED) {
132: if (!da->ghostedcoloring) {
133: PetscMalloc(nc*gnx*gny*sizeof(int),&colors);
134: ii = 0;
135: for (j=gys; j<gys+gny; j++) {
136: for (i=gxs; i<gxs+gnx; 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*gnx*gny,colors,&da->ghostedcoloring);
143: ISColoringSetType(da->ghostedcoloring,IS_COLORING_GHOSTED);
144: }
145: *coloring = da->ghostedcoloring;
146: } else SETERRQ1(1,"Unknown ISColoringType %d",ctype);
147: ISColoringReference(*coloring);
148: return(0);
149: }
151: /* ---------------------------------------------------------------------------------*/
153: int DAGetColoring3d_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring)
154: {
155: int ierr,xs,ys,nx,ny,*colors,i,j,gxs,gys,gnx,gny;
156: int m,n,p,dim,s,k,nc,col,size,zs,gzs,ii,l,nz,gnz;
157: MPI_Comm comm;
158: DAPeriodicType wrap;
159: DAStencilType st;
162: /*
163: nc - number of components per grid point
164: col - number of colors needed in one direction for single component problem
165:
166: */
167: DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);
168: col = 2*s + 1;
169: if (DAXPeriodic(wrap) && (m % col)){
170: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisiblen
171: by 2*stencil_width + 1n");
172: }
173: if (DAYPeriodic(wrap) && (n % col)){
174: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisiblen
175: by 2*stencil_width + 1n");
176: }
177: if (DAZPeriodic(wrap) && (p % col)){
178: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisiblen
179: by 2*stencil_width + 1n");
180: }
182: DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);
183: DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);
184: PetscObjectGetComm((PetscObject)da,&comm);
185: MPI_Comm_size(comm,&size);
187: /* create the coloring */
188: if (ctype == IS_COLORING_LOCAL) {
189: if (!da->localcoloring) {
190: PetscMalloc(nc*nx*ny*nz*sizeof(int),&colors);
191: ii = 0;
192: for (k=zs; k<zs+nz; k++) {
193: for (j=ys; j<ys+ny; j++) {
194: for (i=xs; i<xs+nx; i++) {
195: for (l=0; l<nc; l++) {
196: colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
197: }
198: }
199: }
200: }
201: ISColoringCreate(comm,nc*nx*ny*nz,colors,&da->localcoloring);
202: }
203: *coloring = da->localcoloring;
204: } else if (ctype == IS_COLORING_GHOSTED) {
205: if (!da->ghostedcoloring) {
206: PetscMalloc(nc*gnx*gny*gnz*sizeof(int),&colors);
207: ii = 0;
208: for (k=gzs; k<gzs+gnz; k++) {
209: for (j=gys; j<gys+gny; j++) {
210: for (i=gxs; i<gxs+gnx; i++) {
211: for (l=0; l<nc; l++) {
212: colors[ii++] = l + nc*((i % col) + col*(j % col) + col*col*(k % col));
213: }
214: }
215: }
216: }
217: ISColoringCreate(comm,nc*gnx*gny*gnz,colors,&da->ghostedcoloring);
218: ISColoringSetType(da->ghostedcoloring,IS_COLORING_GHOSTED);
219: }
220: *coloring = da->ghostedcoloring;
221: } else SETERRQ1(1,"Unknown ISColoringType %d",ctype);
222: ISColoringReference(*coloring);
223: return(0);
224: }
226: /* ---------------------------------------------------------------------------------*/
228: int DAGetColoring1d_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring)
229: {
230: int ierr,xs,nx,*colors,i,i1,gxs,gnx,l;
231: int m,dim,s,nc,col,size;
232: MPI_Comm comm;
233: DAPeriodicType wrap;
236: /*
237: nc - number of components per grid point
238: col - number of colors needed in one direction for single component problem
239:
240: */
241: DAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&wrap,0);
242: col = 2*s + 1;
244: if (DAXPeriodic(wrap) && (m % col)) {
245: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points is divisiblen
246: by 2*stencil_width + 1n");
247: }
249: DAGetCorners(da,&xs,0,0,&nx,0,0);
250: DAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);
251: PetscObjectGetComm((PetscObject)da,&comm);
252: MPI_Comm_size(comm,&size);
254: /* create the coloring */
255: if (ctype == IS_COLORING_LOCAL) {
256: if (!da->localcoloring) {
257: PetscMalloc(nc*nx*sizeof(int),&colors);
258: i1 = 0;
259: for (i=xs; i<xs+nx; i++) {
260: for (l=0; l<nc; l++) {
261: colors[i1++] = l + nc*(i % col);
262: }
263: }
264: ISColoringCreate(comm,nc*nx,colors,&da->localcoloring);
265: }
266: *coloring = da->localcoloring;
267: } else if (ctype == IS_COLORING_GHOSTED) {
268: if (!da->ghostedcoloring) {
269: PetscMalloc(nc*gnx*sizeof(int),&colors);
270: i1 = 0;
271: for (i=gxs; i<gxs+gnx; i++) {
272: for (l=0; l<nc; l++) {
273: colors[i1++] = l + nc*(i % col);
274: }
275: }
276: ISColoringCreate(comm,nc*gnx,colors,&da->ghostedcoloring);
277: ISColoringSetType(da->ghostedcoloring,IS_COLORING_GHOSTED);
278: }
279: *coloring = da->ghostedcoloring;
280: } else SETERRQ1(1,"Unknown ISColoringType %d",ctype);
281: ISColoringReference(*coloring);
282: return(0);
283: }
285: int DAGetColoring2d_5pt_MPIAIJ(DA da,ISColoringType ctype,ISColoring *coloring)
286: {
287: int ierr,xs,ys,nx,ny,*colors,i,j,ii,gxs,gys,gnx,gny;
288: int m,n,dim,w,s,k,nc;
289: MPI_Comm comm;
292: /*
293: nc - number of components per grid point
294: col - number of colors needed in one direction for single component problem
295:
296: */
297: ierr = DAGetInfo(da,&dim,&m,&n,0,0,0,0,&w,&s,0,0);
298: nc = w;
299: ierr = DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);
300: ierr = DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);
301: ierr = PetscObjectGetComm((PetscObject)da,&comm);
303: /* create the coloring */
304: if (ctype == IS_COLORING_LOCAL) {
305: if (!da->localcoloring) {
306: PetscMalloc(nc*nx*ny*sizeof(int),&colors);
307: ii = 0;
308: for (j=ys; j<ys+ny; j++) {
309: for (i=xs; i<xs+nx; i++) {
310: for (k=0; k<nc; k++) {
311: colors[ii++] = k + nc*((3*j+i) % 5);
312: }
313: }
314: }
315: ISColoringCreate(comm,nc*nx*ny,colors,&da->localcoloring);
316: }
317: *coloring = da->localcoloring;
318: } else if (ctype == IS_COLORING_GHOSTED) {
319: if (!da->ghostedcoloring) {
320: PetscMalloc(nc*gnx*gny*sizeof(int),&colors);
321: ii = 0;
322: for (j=gys; j<gys+gny; j++) {
323: for (i=gxs; i<gxs+gnx; i++) {
324: for (k=0; k<nc; k++) {
325: colors[ii++] = k + nc*((3*j+i) % 5);
326: }
327: }
328: }
329: ISColoringCreate(comm,nc*gnx*gny,colors,&da->ghostedcoloring);
330: ISColoringSetType(da->ghostedcoloring,IS_COLORING_GHOSTED);
331: }
332: *coloring = da->ghostedcoloring;
333: } else SETERRQ1(1,"Unknown ISColoringType %d",ctype);
334: return(0);
335: }
337: /* =========================================================================== */
338: EXTERN int DAGetMatrix1d_MPIAIJ(DA,Mat *);
339: EXTERN int DAGetMatrix2d_MPIAIJ(DA,Mat *);
340: EXTERN int DAGetMatrix3d_MPIAIJ(DA,Mat *);
341: EXTERN int DAGetMatrix3d_MPIBAIJ(DA,Mat *);
343: /*@C
344: DAGetMatrix - Gets the nonzero structure required for computing the Jacobian via
345: finite differences on a function defined using a stencil on the DA.
347: Collective on DA
349: Input Parameter:
350: + da - the distributed array
351: - mtype - either MATMPIAIJ or MATMPIBAIJ
353: Output Parameters:
354: . J - matrix with the correct nonzero structure
355: (obviously without the correct Jacobian values)
357: Level: advanced
360: This does not yet handle BAIJ matrices, because
361: 1) we need a way for the user to indicate what matrix type they want
362: 2) after the matrix is created, for BAIJ matrices we need to set nc to 1 and
363: use MatSetValuesBlockedLocal() instead of MatSetValuesLocal()
365: .seealso ISColoringView(), ISColoringGetIS(), MatFDColoringCreate(), DAGetMatrixMPIBAIJ()
367: @*/
368: int DAGetMatrix(DA da,MatType mtype,Mat *J)
369: {
370: int ierr,dim;
371: PetscTruth aij;
374: /*
375: m
376: ------------------------------------------------------
377: | |
378: | |
379: | ---------------------- |
380: | | | |
381: n | yn | | |
382: | | | |
383: | .--------------------- |
384: | (xs,ys) xn |
385: | . |
386: | (gxs,gys) |
387: | |
388: -----------------------------------------------------
389: */
391: /*
392: nc - number of components per grid point
393: col - number of colors needed in one direction for single component problem
394:
395: */
396: DAGetInfo(da,&dim,0,0,0,0,0,0,0,0,0,0);
398: /*
399: We do not provide a getcoloring function in the DA operations because
400: the basic DA does not know about matrices. We think of DA as being more
401: more low-level then matrices.
402: */
403: PetscStrcmp(MATMPIAIJ,mtype,&aij);
404: if (aij) {
405: if (dim == 1) {
406: DAGetMatrix1d_MPIAIJ(da,J);
407: } else if (dim == 2) {
408: DAGetMatrix2d_MPIAIJ(da,J);
409: } else if (dim == 3) {
410: DAGetMatrix3d_MPIAIJ(da,J);
411: }
412: } else {
413: if (dim == 3) {
414: DAGetMatrix3d_MPIBAIJ(da,J);
415: } else {
416: SETERRQ1(1,"Not done for %d dimension, send us mail petsc-maint@mcs.anl.gov for code",dim);
417: }
418: }
419: return(0);
420: }
422: /* ---------------------------------------------------------------------------------*/
424: int DAGetMatrix2d_MPIAIJ(DA da,Mat *J)
425: {
426: int ierr,xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
427: int m,n,dim,s,*cols,k,nc,*rows,col,cnt,l,p;
428: int lstart,lend,pstart,pend,*dnz,*onz,size;
429: int dims[2],starts[2];
430: MPI_Comm comm;
431: PetscScalar *values;
432: DAPeriodicType wrap;
433: ISLocalToGlobalMapping ltog;
434: DAStencilType st;
437: /*
438: nc - number of components per grid point
439: col - number of colors needed in one direction for single component problem
440:
441: */
442: DAGetInfo(da,&dim,&m,&n,0,0,0,0,&nc,&s,&wrap,&st);
443: col = 2*s + 1;
444: if (DAXPeriodic(wrap) && (m % col)){
445: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisiblen
446: by 2*stencil_width + 1n");
447: }
448: if (DAYPeriodic(wrap) && (n % col)){
449: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisiblen
450: by 2*stencil_width + 1n");
451: }
452: DAGetCorners(da,&xs,&ys,0,&nx,&ny,0);
453: DAGetGhostCorners(da,&gxs,&gys,0,&gnx,&gny,0);
454: PetscObjectGetComm((PetscObject)da,&comm);
455: MPI_Comm_size(comm,&size);
457: /* create empty Jacobian matrix */
458: ierr = MatCreate(comm,nc*nx*ny,nc*nx*ny,PETSC_DECIDE,PETSC_DECIDE,J);
460: PetscMalloc(col*col*nc*nc*sizeof(PetscScalar),&values);
461: PetscMemzero(values,col*col*nc*nc*sizeof(PetscScalar));
462: PetscMalloc(nc*sizeof(int),&rows);
463: PetscMalloc(col*col*nc*nc*sizeof(int),&cols);
464: DAGetISLocalToGlobalMapping(da,<og);
465:
466: /* determine the matrix preallocation information */
467: MatPreallocateInitialize(comm,nc*nx*ny,nc*nx*ny,dnz,onz);
468: for (i=xs; i<xs+nx; i++) {
470: pstart = PetscMax(-s,-i);
471: pend = PetscMin(s,m-i-1);
473: for (j=ys; j<ys+ny; j++) {
474: slot = i - gxs + gnx*(j - gys);
476: lstart = PetscMax(-s,-j);
477: lend = PetscMin(s,n-j-1);
479: cnt = 0;
480: for (k=0; k<nc; k++) {
481: for (l=lstart; l<lend+1; l++) {
482: for (p=pstart; p<pend+1; p++) {
483: if ((st == DA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */
484: cols[cnt++] = k + nc*(slot + gnx*l + p);
485: }
486: }
487: }
488: rows[k] = k + nc*(slot);
489: }
490: MatPreallocateSetLocal(ltog,nc,rows,cnt,cols,dnz,onz);
491: }
492: }
493: /* set matrix type and preallocation information */
494: if (size > 1) {
495: MatSetType(*J,MATMPIAIJ);
496: } else {
497: MatSetType(*J,MATSEQAIJ);
498: }
499: MatSeqAIJSetPreallocation(*J,0,dnz);
500: MatSeqBAIJSetPreallocation(*J,nc,0,dnz);
501: MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);
502: MatMPIBAIJSetPreallocation(*J,nc,0,dnz,0,onz);
503: MatPreallocateFinalize(dnz,onz);
504: MatSetLocalToGlobalMapping(*J,ltog);
505: DAGetGhostCorners(da,&starts[0],&starts[1],PETSC_IGNORE,&dims[0],&dims[1],PETSC_IGNORE);
506: MatSetStencil(*J,2,dims,starts,nc);
508: /*
509: For each node in the grid: we get the neighbors in the local (on processor ordering
510: that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
511: PETSc ordering.
512: */
513: for (i=xs; i<xs+nx; i++) {
514:
515: pstart = PetscMax(-s,-i);
516: pend = PetscMin(s,m-i-1);
517:
518: for (j=ys; j<ys+ny; j++) {
519: slot = i - gxs + gnx*(j - gys);
520:
521: lstart = PetscMax(-s,-j);
522: lend = PetscMin(s,n-j-1);
524: cnt = 0;
525: for (k=0; k<nc; k++) {
526: for (l=lstart; l<lend+1; l++) {
527: for (p=pstart; p<pend+1; p++) {
528: if ((st == DA_STENCIL_BOX) || (!l || !p)) { /* entries on star have either l = 0 or p = 0 */
529: cols[cnt++] = k + nc*(slot + gnx*l + p);
530: }
531: }
532: }
533: rows[k] = k + nc*(slot);
534: }
535: MatSetValuesLocal(*J,nc,rows,cnt,cols,values,INSERT_VALUES);
536: }
537: }
538: PetscFree(values);
539: PetscFree(rows);
540: PetscFree(cols);
541: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
542: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
543: return(0);
544: }
546: /* ---------------------------------------------------------------------------------*/
548: int DAGetMatrix3d_MPIAIJ(DA da,Mat *J)
549: {
550: int ierr,xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
551: int m,n,dim,s,*cols,k,nc,*rows,col,cnt,l,p,*dnz,*onz;
552: int istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk,size;
553: int dims[3],starts[3];
554: MPI_Comm comm;
555: PetscScalar *values;
556: DAPeriodicType wrap;
557: ISLocalToGlobalMapping ltog;
558: DAStencilType st;
561: /*
562: nc - number of components per grid point
563: col - number of colors needed in one direction for single component problem
564:
565: */
566: DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);
567: col = 2*s + 1;
568: if (DAXPeriodic(wrap) && (m % col)){
569: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in X is divisiblen
570: by 2*stencil_width + 1n");
571: }
572: if (DAYPeriodic(wrap) && (n % col)){
573: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Y is divisiblen
574: by 2*stencil_width + 1n");
575: }
576: if (DAZPeriodic(wrap) && (p % col)){
577: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points in Z is divisiblen
578: by 2*stencil_width + 1n");
579: }
581: DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);
582: DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);
583: PetscObjectGetComm((PetscObject)da,&comm);
584: MPI_Comm_size(comm,&size);
587: /* create the matrix */
588: /* create empty Jacobian matrix */
589: MatCreate(comm,nc*nx*ny*nz,nc*nx*ny*nz,PETSC_DECIDE,PETSC_DECIDE,J);
590: PetscMalloc(col*col*col*nc*nc*nc*sizeof(PetscScalar),&values);
591: PetscMemzero(values,col*col*col*nc*nc*nc*sizeof(PetscScalar));
592: PetscMalloc(nc*sizeof(int),&rows);
593: PetscMalloc(col*col*col*nc*sizeof(int),&cols);
594: DAGetISLocalToGlobalMapping(da,<og);
596: /* determine the matrix preallocation information */
597: MatPreallocateInitialize(comm,nc*nx*ny*nz,nc*nx*ny*nz,dnz,onz);
598: for (i=xs; i<xs+nx; i++) {
599: istart = PetscMax(-s,-i);
600: iend = PetscMin(s,m-i-1);
601: for (j=ys; j<ys+ny; j++) {
602: jstart = PetscMax(-s,-j);
603: jend = PetscMin(s,n-j-1);
604: for (k=zs; k<zs+nz; k++) {
605: kstart = PetscMax(-s,-k);
606: kend = PetscMin(s,p-k-1);
607:
608: slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
609:
610: cnt = 0;
611: for (l=0; l<nc; l++) {
612: for (ii=istart; ii<iend+1; ii++) {
613: for (jj=jstart; jj<jend+1; jj++) {
614: for (kk=kstart; kk<kend+1; kk++) {
615: if ((st == DA_STENCIL_BOX) || (!ii || !jj || !kk)) { /* entries on star */
616: cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
617: }
618: }
619: }
620: }
621: rows[l] = l + nc*(slot);
622: }
623: MatPreallocateSetLocal(ltog,nc,rows,cnt,cols,dnz,onz);
624: }
625: }
626: }
627: /* set matrix type and preallocation */
628: if (size > 1) {
629: MatSetType(*J,MATMPIAIJ);
630: } else {
631: MatSetType(*J,MATSEQAIJ);
632: }
633: MatSeqAIJSetPreallocation(*J,0,dnz);
634: MatSeqBAIJSetPreallocation(*J,nc,0,dnz);
635: MatMPIAIJSetPreallocation(*J,0,dnz,0,onz);
636: MatMPIBAIJSetPreallocation(*J,nc,0,dnz,0,onz);
637: MatPreallocateFinalize(dnz,onz);
638: MatSetLocalToGlobalMapping(*J,ltog);
639: DAGetGhostCorners(da,&starts[0],&starts[1],&starts[2],&dims[0],&dims[1],&dims[2]);
640: MatSetStencil(*J,3,dims,starts,nc);
642: /*
643: For each node in the grid: we get the neighbors in the local (on processor ordering
644: that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
645: PETSc ordering.
646: */
647: for (i=xs; i<xs+nx; i++) {
648: istart = PetscMax(-s,-i);
649: iend = PetscMin(s,m-i-1);
650: for (j=ys; j<ys+ny; j++) {
651: jstart = PetscMax(-s,-j);
652: jend = PetscMin(s,n-j-1);
653: for (k=zs; k<zs+nz; k++) {
654: kstart = PetscMax(-s,-k);
655: kend = PetscMin(s,p-k-1);
656:
657: slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
658:
659: cnt = 0;
660: for (l=0; l<nc; l++) {
661: for (ii=istart; ii<iend+1; ii++) {
662: for (jj=jstart; jj<jend+1; jj++) {
663: for (kk=kstart; kk<kend+1; kk++) {
664: if ((st == DA_STENCIL_BOX) || (!ii || !jj || !kk)) { /* entries on star */
665: cols[cnt++] = l + nc*(slot + ii + gnx*jj + gnx*gny*kk);
666: }
667: }
668: }
669: }
670: rows[l] = l + nc*(slot);
671: }
672: MatSetValuesLocal(*J,nc,rows,cnt,cols,values,INSERT_VALUES);
673: }
674: }
675: }
676: PetscFree(values);
677: PetscFree(rows);
678: PetscFree(cols);
679: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
680: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
681: return(0);
682: }
684: /* ---------------------------------------------------------------------------------*/
686: int DAGetMatrix1d_MPIAIJ(DA da,Mat *J)
687: {
688: int ierr,xs,nx,i,i1,slot,gxs,gnx;
689: int m,dim,s,*cols,nc,*rows,col,cnt,l;
690: int istart,iend,size;
691: int dims[1],starts[1];
692: MPI_Comm comm;
693: PetscScalar *values;
694: DAPeriodicType wrap;
695: ISLocalToGlobalMapping ltog;
698: /*
699: nc - number of components per grid point
700: col - number of colors needed in one direction for single component problem
701:
702: */
703: DAGetInfo(da,&dim,&m,0,0,0,0,0,&nc,&s,&wrap,0);
704: col = 2*s + 1;
706: if (DAXPeriodic(wrap) && (m % col)) {
707: SETERRQ(PETSC_ERR_SUP,"For coloring efficiency ensure number of grid points is divisiblen
708: by 2*stencil_width + 1n");
709: }
712: DAGetCorners(da,&xs,0,0,&nx,0,0);
713: DAGetGhostCorners(da,&gxs,0,0,&gnx,0,0);
714: PetscObjectGetComm((PetscObject)da,&comm);
715: MPI_Comm_size(comm,&size);
717: /* create empty Jacobian matrix */
718:
719: ierr = MatCreate(comm,nc*nx,nc*nx,PETSC_DECIDE,PETSC_DECIDE,J);
720: if (size > 1) {
721: MatSetType(*J,MATMPIAIJ);
722: } else {
723: MatSetType(*J,MATSEQAIJ);
724: }
725: MatSeqAIJSetPreallocation(*J,col*nc,0);
726: MatSeqBAIJSetPreallocation(*J,nc,col,0);
727: MatMPIAIJSetPreallocation(*J,col*nc,0,0,0);
728: MatMPIBAIJSetPreallocation(*J,nc,col,0,0,0);
729: DAGetGhostCorners(da,&starts[0],PETSC_IGNORE,PETSC_IGNORE,&dims[0],PETSC_IGNORE,PETSC_IGNORE);
730: MatSetStencil(*J,1,dims,starts,nc);
731:
732: PetscMalloc(col*nc*nc*sizeof(PetscScalar),&values);
733: PetscMemzero(values,col*nc*nc*sizeof(PetscScalar));
734: PetscMalloc(nc*sizeof(int),&rows);
735: PetscMalloc(col*nc*sizeof(int),&cols);
736:
737: DAGetISLocalToGlobalMapping(da,<og);
738: MatSetLocalToGlobalMapping(*J,ltog);
739:
740: /*
741: For each node in the grid: we get the neighbors in the local (on processor ordering
742: that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
743: PETSc ordering.
744: */
745: for (i=xs; i<xs+nx; i++) {
746: istart = PetscMax(-s,gxs - i);
747: iend = PetscMin(s,gxs + gnx - i - 1);
748: slot = i - gxs;
749:
750: cnt = 0;
751: for (l=0; l<nc; l++) {
752: for (i1=istart; i1<iend+1; i1++) {
753: cols[cnt++] = l + nc*(slot + i1);
754: }
755: rows[l] = l + nc*(slot);
756: }
757: MatSetValuesLocal(*J,nc,rows,cnt,cols,values,INSERT_VALUES);
758: }
759: PetscFree(values);
760: PetscFree(rows);
761: PetscFree(cols);
762: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
763: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
764: return(0);
765: }
767: int DAGetMatrix3d_MPIBAIJ(DA da,Mat *J)
768: {
769: int ierr,xs,ys,nx,ny,i,j,slot,gxs,gys,gnx,gny;
770: int m,n,dim,s,*cols,k,nc,col,cnt,p,*dnz,*onz;
771: int istart,iend,jstart,jend,kstart,kend,zs,nz,gzs,gnz,ii,jj,kk;
772: MPI_Comm comm;
773: PetscScalar *values;
774: DAPeriodicType wrap;
775: DAStencilType st;
776: ISLocalToGlobalMapping ltog;
779: /*
780: nc - number of components per grid point
781: col - number of colors needed in one direction for single component problem
782:
783: */
784: DAGetInfo(da,&dim,&m,&n,&p,0,0,0,&nc,&s,&wrap,&st);
785: if (wrap != DA_NONPERIODIC) SETERRQ(PETSC_ERR_SUP,"Currently no support for periodic");
786: col = 2*s + 1;
788: DAGetCorners(da,&xs,&ys,&zs,&nx,&ny,&nz);
789: DAGetGhostCorners(da,&gxs,&gys,&gzs,&gnx,&gny,&gnz);
790: PetscObjectGetComm((PetscObject)da,&comm);
792: /* create the matrix */
793: ierr = PetscMalloc(col*col*col*nc*nc*sizeof(PetscScalar),&values);
794: ierr = PetscMemzero(values,col*col*col*nc*nc*sizeof(PetscScalar));
795: ierr = PetscMalloc(col*col*col*sizeof(int),&cols);
797: DAGetISLocalToGlobalMappingBlck(da,<og);
799: /* determine the matrix preallocation information */
800: MatPreallocateInitialize(comm,nx*ny*nz,nx*ny*nz,dnz,onz);
801: for (i=xs; i<xs+nx; i++) {
802: istart = PetscMax(-s,-i);
803: iend = PetscMin(s,m-i-1);
804: for (j=ys; j<ys+ny; j++) {
805: jstart = PetscMax(-s,-j);
806: jend = PetscMin(s,n-j-1);
807: for (k=zs; k<zs+nz; k++) {
808: kstart = PetscMax(-s,-k);
809: kend = PetscMin(s,p-k-1);
811: slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
813: /* Find block columns in block row */
814: cnt = 0;
815: if (st == DA_STENCIL_BOX) { /* if using BOX stencil */
816: for (ii=istart; ii<iend+1; ii++) {
817: for (jj=jstart; jj<jend+1; jj++) {
818: for (kk=kstart; kk<kend+1; kk++) {
819: cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
820: }
821: }
822: }
823: } else { /* Star stencil */
824: cnt = 0;
825: for (ii=istart; ii<iend+1; ii++) {
826: if (ii) {
827: /* jj and kk must be zero */
828: /* cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; */
829: cols[cnt++] = slot + ii;
830: } else {
831: for (jj=jstart; jj<jend+1; jj++) {
832: if (jj) {
833: /* ii and kk must be zero */
834: cols[cnt++] = slot + gnx*jj;
835: } else {
836: /* ii and jj must be zero */
837: for (kk=kstart; kk<kend+1; kk++) {
838: cols[cnt++] = slot + gnx*gny*kk;
839: }
840: }
841: }
842: }
843: }
844: }
845: MatPreallocateSetLocal(ltog,1,&slot,cnt,cols,dnz,onz);
846: }
847: }
848: }
850: /* create empty Jacobian matrix */
851: MatCreateMPIBAIJ(comm,nc,nc*nx*ny*nz,nc*nx*ny*nz,PETSC_DECIDE,
852: PETSC_DECIDE,0,dnz,0,onz,J);
854: MatPreallocateFinalize(dnz,onz);
855: MatSetLocalToGlobalMappingBlock(*J,ltog);
857: /*
858: For each node in the grid: we get the neighbors in the local (on processor ordering
859: that includes the ghost points) then MatSetValuesLocal() maps those indices to the global
860: PETSc ordering.
861: */
863: for (i=xs; i<xs+nx; i++) {
864: istart = PetscMax(-s,-i);
865: iend = PetscMin(s,m-i-1);
866: for (j=ys; j<ys+ny; j++) {
867: jstart = PetscMax(-s,-j);
868: jend = PetscMin(s,n-j-1);
869: for (k=zs; k<zs+nz; k++) {
870: kstart = PetscMax(-s,-k);
871: kend = PetscMin(s,p-k-1);
872:
873: slot = i - gxs + gnx*(j - gys) + gnx*gny*(k - gzs);
874:
875: cnt = 0;
876: if (st == DA_STENCIL_BOX) { /* if using BOX stencil */
877: for (ii=istart; ii<iend+1; ii++) {
878: for (jj=jstart; jj<jend+1; jj++) {
879: for (kk=kstart; kk<kend+1; kk++) {
880: cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk;
881: }
882: }
883: }
884: } else { /* Star stencil */
885: cnt = 0;
886: for (ii=istart; ii<iend+1; ii++) {
887: if (ii) {
888: /* jj and kk must be zero */
889: /* cols[cnt++] = slot + ii + gnx*jj + gnx*gny*kk; */
890: cols[cnt++] = slot + ii;
891: } else {
892: for (jj=jstart; jj<jend+1; jj++) {
893: if (jj) {
894: /* ii and kk must be zero */
895: cols[cnt++] = slot + gnx*jj;
896: } else {
897: /* ii and jj must be zero */
898: for (kk=kstart; kk<kend+1; kk++) {
899: cols[cnt++] = slot + gnx*gny*kk;
900: }
901: }
902: }
903: }
904: }
905: }
906: MatSetValuesBlockedLocal(*J,1,&slot,cnt,cols,values,INSERT_VALUES);
907: }
908: }
909: }
910: PetscFree(values);
911: PetscFree(cols);
912: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
913: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
914: return(0);
915: }