Actual source code: da3.c
1: /*$Id: da3.c,v 1.136 2001/09/07 20:12:17 bsmith Exp $*/
3: /*
4: Code for manipulating distributed regular 3d arrays in parallel.
5: File created by Peter Mell 7/14/95
6: */
8: #include src/dm/da/daimpl.h
10: #if defined (PETSC_HAVE_AMS)
11: EXTERN_C_BEGIN
12: EXTERN int AMSSetFieldBlock_DA(AMS_Memory,char *,Vec);
13: EXTERN_C_END
14: #endif
18: int DAView_3d(DA da,PetscViewer viewer)
19: {
20: int rank,ierr;
21: PetscTruth isascii,isdraw;
24: MPI_Comm_rank(da->comm,&rank);
26: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
27: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_DRAW,&isdraw);
28: if (isascii) {
29: PetscViewerASCIISynchronizedPrintf(viewer,"Processor [%d] M %d N %d P %d m %d n %d p %d w %d s %d\n",
30: rank,da->M,da->N,da->P,da->m,da->n,da->p,da->w,da->s);
31: PetscViewerASCIISynchronizedPrintf(viewer,"X range of indices: %d %d, Y range of indices: %d %d, Z range of indices: %d %d\n",
32: da->xs,da->xe,da->ys,da->ye,da->zs,da->ze);
33: #if !defined(PETSC_USE_COMPLEX)
34: if (da->coordinates) {
35: int last;
36: PetscReal *coors;
37: VecGetArray(da->coordinates,&coors);
38: VecGetLocalSize(da->coordinates,&last);
39: last = last - 3;
40: PetscViewerASCIISynchronizedPrintf(viewer,"Lower left corner %g %g %g : Upper right %g %g %g\n",
41: coors[0],coors[1],coors[2],coors[last],coors[last+1],coors[last+2]);
42: VecRestoreArray(da->coordinates,&coors);
43: }
44: #endif
45: PetscViewerFlush(viewer);
46: } else if (isdraw) {
47: PetscDraw draw;
48: PetscReal ymin = -1.0,ymax = (PetscReal)da->N;
49: PetscReal xmin = -1.0,xmax = (PetscReal)((da->M+2)*da->P),x,y,ycoord,xcoord;
50: int k,plane,base,*idx;
51: char node[10];
52: PetscTruth isnull;
54: PetscViewerDrawGetDraw(viewer,0,&draw);
55: PetscDrawIsNull(draw,&isnull); if (isnull) return(0);
56: PetscDrawSetCoordinates(draw,xmin,ymin,xmax,ymax);
57: PetscDrawSynchronizedClear(draw);
59: /* first processor draw all node lines */
60: if (!rank) {
61: for (k=0; k<da->P; k++) {
62: ymin = 0.0; ymax = (PetscReal)(da->N - 1);
63: for (xmin=(PetscReal)(k*(da->M+1)); xmin<(PetscReal)(da->M+(k*(da->M+1))); xmin++) {
64: PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_BLACK);
65: }
66:
67: xmin = (PetscReal)(k*(da->M+1)); xmax = xmin + (PetscReal)(da->M - 1);
68: for (ymin=0; ymin<(PetscReal)da->N; ymin++) {
69: PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_BLACK);
70: }
71: }
72: }
73: PetscDrawSynchronizedFlush(draw);
74: PetscDrawPause(draw);
76: for (k=0; k<da->P; k++) { /*Go through and draw for each plane*/
77: if ((k >= da->zs) && (k < da->ze)) {
78: /* draw my box */
79: ymin = da->ys;
80: ymax = da->ye - 1;
81: xmin = da->xs/da->w + (da->M+1)*k;
82: xmax =(da->xe-1)/da->w + (da->M+1)*k;
84: PetscDrawLine(draw,xmin,ymin,xmax,ymin,PETSC_DRAW_RED);
85: PetscDrawLine(draw,xmin,ymin,xmin,ymax,PETSC_DRAW_RED);
86: PetscDrawLine(draw,xmin,ymax,xmax,ymax,PETSC_DRAW_RED);
87: PetscDrawLine(draw,xmax,ymin,xmax,ymax,PETSC_DRAW_RED);
89: xmin = da->xs/da->w;
90: xmax =(da->xe-1)/da->w;
92: /* put in numbers*/
93: base = (da->base+(da->xe-da->xs)*(da->ye-da->ys)*(k-da->zs))/da->w;
95: /* Identify which processor owns the box */
96: sprintf(node,"%d",rank);
97: PetscDrawString(draw,xmin+(da->M+1)*k+.2,ymin+.3,PETSC_DRAW_RED,node);
99: for (y=ymin; y<=ymax; y++) {
100: for (x=xmin+(da->M+1)*k; x<=xmax+(da->M+1)*k; x++) {
101: sprintf(node,"%d",base++);
102: PetscDrawString(draw,x,y,PETSC_DRAW_BLACK,node);
103: }
104: }
105:
106: }
107: }
108: PetscDrawSynchronizedFlush(draw);
109: PetscDrawPause(draw);
111: for (k=0-da->s; k<da->P+da->s; k++) {
112: /* Go through and draw for each plane */
113: if ((k >= da->Zs) && (k < da->Ze)) {
114:
115: /* overlay ghost numbers, useful for error checking */
116: base = (da->Xe-da->Xs)*(da->Ye-da->Ys)*(k-da->Zs); idx = da->idx;
117: plane=k;
118: /* Keep z wrap around points on the dradrawg */
119: if (k<0) { plane=da->P+k; }
120: if (k>=da->P) { plane=k-da->P; }
121: ymin = da->Ys; ymax = da->Ye;
122: xmin = (da->M+1)*plane*da->w;
123: xmax = (da->M+1)*plane*da->w+da->M*da->w;
124: for (y=ymin; y<ymax; y++) {
125: for (x=xmin+da->Xs; x<xmin+da->Xe; x+=da->w) {
126: sprintf(node,"%d",idx[base]/da->w);
127: ycoord = y;
128: /*Keep y wrap around points on drawing */
129: if (y<0) { ycoord = da->N+y; }
131: if (y>=da->N) { ycoord = y-da->N; }
132: xcoord = x; /* Keep x wrap points on drawing */
134: if (x<xmin) { xcoord = xmax - (xmin-x); }
135: if (x>=xmax) { xcoord = xmin + (x-xmax); }
136: PetscDrawString(draw,xcoord/da->w,ycoord,PETSC_DRAW_BLUE,node);
137: base+=da->w;
138: }
139: }
140: }
141: }
142: PetscDrawSynchronizedFlush(draw);
143: PetscDrawPause(draw);
144: } else {
145: SETERRQ1(1,"Viewer type %s not supported for DA 3d",((PetscObject)viewer)->type_name);
146: }
147: return(0);
148: }
150: EXTERN int DAPublish_Petsc(PetscObject);
154: /*@C
155: DACreate3d - Creates an object that will manage the communication of three-dimensional
156: regular array data that is distributed across some processors.
158: Collective on MPI_Comm
160: Input Parameters:
161: + comm - MPI communicator
162: . wrap - type of periodicity the array should have, if any. Use one
163: of DA_NONPERIODIC, DA_XPERIODIC, DA_YPERIODIC, DA_XYPERIODIC, DA_XYZPERIODIC, DA_XZPERIODIC, or DA_YZPERIODIC.
164: . stencil_type - Type of stencil (DA_STENCIL_STAR or DA_STENCIL_BOX)
165: . M,N,P - global dimension in each direction of the array (use -M, -N, and or -P to indicate that it may be set to a different value
166: from the command line with -da_grid_x <M> -da_grid_y <N> -da_grid_z <P>)
167: . m,n,p - corresponding number of processors in each dimension
168: (or PETSC_DECIDE to have calculated)
169: . dof - number of degrees of freedom per node
170: . lx, ly, lz - arrays containing the number of nodes in each cell along
171: the x, y, and z coordinates, or PETSC_NULL. If non-null, these
172: must be of length as m,n,p and the corresponding
173: m,n, or p cannot be PETSC_DECIDE. Sum of the lx[] entries must be M, sum of
174: the ly[] must n, sum of the lz[] must be P
175: - s - stencil width
177: Output Parameter:
178: . inra - the resulting distributed array object
180: Options Database Key:
181: + -da_view - Calls DAView() at the conclusion of DACreate3d()
182: . -da_grid_x <nx> - number of grid points in x direction, if M < 0
183: . -da_grid_y <ny> - number of grid points in y direction, if N < 0
184: - -da_grid_z <nz> - number of grid points in z direction, if P < 0
186: Level: beginner
188: Notes:
189: The stencil type DA_STENCIL_STAR with width 1 corresponds to the
190: standard 7-pt stencil, while DA_STENCIL_BOX with width 1 denotes
191: the standard 27-pt stencil.
193: The array data itself is NOT stored in the DA, it is stored in Vec objects;
194: The appropriate vector objects can be obtained with calls to DACreateGlobalVector()
195: and DACreateLocalVector() and calls to VecDuplicate() if more are needed.
197: .keywords: distributed array, create, three-dimensional
199: .seealso: DADestroy(), DAView(), DACreate1d(), DACreate2d(), DAGlobalToLocalBegin(),
200: DAGlobalToLocalEnd(), DALocalToGlobal(), DALocalToLocalBegin(), DALocalToLocalEnd(),
201: DAGetInfo(), DACreateGlobalVector(), DACreateLocalVector(), DACreateNaturalVector(), DALoad(), DAView()
203: @*/
204: int DACreate3d(MPI_Comm comm,DAPeriodicType wrap,DAStencilType stencil_type,int M,
205: int N,int P,int m,int n,int p,int dof,int s,int *lx,int *ly,int *lz,DA *inra)
206: {
207: int rank,size,ierr,start,end,pm;
208: int xs,xe,ys,ye,zs,ze,x,y,z,Xs,Xe,Ys,Ye,Zs,Ze;
209: int left,up,down,bottom,top,i,j,k,*idx,nn,*flx = 0,*fly = 0,*flz = 0;
210: int n0,n1,n2,n3,n4,n5,n6,n7,n8,n9,n10,n11,n12,n14;
211: int n15,n16,n17,n18,n19,n20,n21,n22,n23,n24,n25,n26;
212: int *bases,*ldims,x_t,y_t,z_t,s_t,base,count,s_x,s_y,s_z;
213: int tM = M,tN = N,tP = P;
214: int sn0 = 0,sn1 = 0,sn2 = 0,sn3 = 0,sn5 = 0,sn6 = 0,sn7 = 0;
215: int sn8 = 0,sn9 = 0,sn11 = 0,sn15 = 0,sn24 = 0,sn25 = 0,sn26 = 0;
216: int sn17 = 0,sn18 = 0,sn19 = 0,sn20 = 0,sn21 = 0,sn23 = 0,refine_x = 2, refine_y = 2, refine_z = 2;
217: PetscTruth flg1,flg2;
218: DA da;
219: Vec local,global;
220: VecScatter ltog,gtol;
221: IS to,from;
225: *inra = 0;
226: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
227: DMInitializePackage(PETSC_NULL);
228: #endif
230: if (dof < 1) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Must have 1 or more degrees of freedom per node: %d",dof);
231: if (s < 0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Stencil width cannot be negative: %d",s);
233: PetscOptionsBegin(comm,PETSC_NULL,"3d DA Options","DA");
234: if (M < 0){
235: tM = -M;
236: PetscOptionsInt("-da_grid_x","Number of grid points in x direction","DACreate3d",tM,&tM,PETSC_NULL);
237: }
238: if (N < 0){
239: tN = -N;
240: PetscOptionsInt("-da_grid_y","Number of grid points in y direction","DACreate3d",tN,&tN,PETSC_NULL);
241: }
242: if (P < 0){
243: tP = -P;
244: PetscOptionsInt("-da_grid_z","Number of grid points in z direction","DACreate3d",tP,&tP,PETSC_NULL);
245: }
246: PetscOptionsInt("-da_processors_x","Number of processors in x direction","DACreate3d",m,&m,PETSC_NULL);
247: PetscOptionsInt("-da_processors_y","Number of processors in y direction","DACreate3d",n,&n,PETSC_NULL);
248: PetscOptionsInt("-da_processors_z","Number of processors in z direction","DACreate3d",p,&p,PETSC_NULL);
249: PetscOptionsInt("-da_refine_x","Refinement ratio in x direction","DACreate3d",refine_x,&refine_x,PETSC_NULL);
250: PetscOptionsInt("-da_refine_y","Refinement ratio in y direction","DACreate3d",refine_y,&refine_y,PETSC_NULL);
251: PetscOptionsInt("-da_refine_z","Refinement ratio in z direction","DACreate3d",refine_z,&refine_z,PETSC_NULL);
252: PetscOptionsEnd();
253: M = tM; N = tN; P = tP;
255: PetscHeaderCreate(da,_p_DA,struct _DAOps,DA_COOKIE,0,"DA",comm,DADestroy,DAView);
256: da->bops->publish = DAPublish_Petsc;
257: da->ops->createglobalvector = DACreateGlobalVector;
258: da->ops->getinterpolation = DAGetInterpolation;
259: da->ops->getcoloring = DAGetColoring;
260: da->ops->getmatrix = DAGetMatrix;
261: da->ops->refine = DARefine;
263: PetscLogObjectCreate(da);
264: PetscLogObjectMemory(da,sizeof(struct _p_DA));
265: da->dim = 3;
266: da->interptype = DA_Q1;
267: da->refine_x = refine_x;
268: da->refine_y = refine_y;
269: da->refine_z = refine_z;
270: PetscMalloc(dof*sizeof(char*),&da->fieldname);
271: PetscMemzero(da->fieldname,dof*sizeof(char*));
273: MPI_Comm_size(comm,&size);
274: MPI_Comm_rank(comm,&rank);
276: if (m != PETSC_DECIDE) {
277: if (m < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in X direction: %d",m);}
278: else if (m > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in X direction: %d %d",m,size);}
279: }
280: if (n != PETSC_DECIDE) {
281: if (n < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Y direction: %d",n);}
282: else if (n > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Y direction: %d %d",n,size);}
283: }
284: if (p != PETSC_DECIDE) {
285: if (p < 1) {SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Non-positive number of processors in Z direction: %d",p);}
286: else if (p > size) {SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Too many processors in Z direction: %d %d",p,size);}
287: }
289: /* Partition the array among the processors */
290: if (m == PETSC_DECIDE && n != PETSC_DECIDE && p != PETSC_DECIDE) {
291: m = size/(n*p);
292: } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
293: n = size/(m*p);
294: } else if (m != PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
295: p = size/(m*n);
296: } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p != PETSC_DECIDE) {
297: /* try for squarish distribution */
298: m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)N*p)));
299: if (!m) m = 1;
300: while (m > 0) {
301: n = size/(m*p);
302: if (m*n*p == size) break;
303: m--;
304: }
305: if (!m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad p value: p = %d",p);
306: if (M > N && m < n) {int _m = m; m = n; n = _m;}
307: } else if (m == PETSC_DECIDE && n != PETSC_DECIDE && p == PETSC_DECIDE) {
308: /* try for squarish distribution */
309: m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
310: if (!m) m = 1;
311: while (m > 0) {
312: p = size/(m*n);
313: if (m*n*p == size) break;
314: m--;
315: }
316: if (!m) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad n value: n = %d",n);
317: if (M > P && m < p) {int _m = m; m = p; p = _m;}
318: } else if (m != PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
319: /* try for squarish distribution */
320: n = (int)(0.5 + sqrt(((PetscReal)N)*((PetscReal)size)/((PetscReal)P*m)));
321: if (!n) n = 1;
322: while (n > 0) {
323: p = size/(m*n);
324: if (m*n*p == size) break;
325: n--;
326: }
327: if (!n) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"bad m value: m = %d",n);
328: if (N > P && n < p) {int _n = n; n = p; p = _n;}
329: } else if (m == PETSC_DECIDE && n == PETSC_DECIDE && p == PETSC_DECIDE) {
330: /* try for squarish distribution */
331: n = (int)(0.5 + pow(((PetscReal)N*N)*((PetscReal)size)/((PetscReal)P*M),1./3.));
332: if (!n) n = 1;
333: while (n > 0) {
334: pm = size/n;
335: if (n*pm == size) break;
336: n--;
337: }
338: if (!n) n = 1;
339: m = (int)(0.5 + sqrt(((PetscReal)M)*((PetscReal)size)/((PetscReal)P*n)));
340: if (!m) m = 1;
341: while (m > 0) {
342: p = size/(m*n);
343: if (m*n*p == size) break;
344: m--;
345: }
346: if (M > P && m < p) {int _m = m; m = p; p = _m;}
347: } else if (m*n*p != size) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Given Bad partition");
349: if (m*n*p != size) SETERRQ(PETSC_ERR_PLIB,"Could not find good partition");
350: if (M < m) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in x direction is too fine! %d %d",M,m);
351: if (N < n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in y direction is too fine! %d %d",N,n);
352: if (P < p) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Partition in z direction is too fine! %d %d",P,p);
354: PetscOptionsHasName(PETSC_NULL,"-da_partition_nodes_at_end",&flg2);
355: /*
356: Determine locally owned region
357: [x, y, or z]s is the first local node number, [x, y, z] is the number of local nodes
358: */
359: if (lx) { /* user decided distribution */
360: x = lx[rank % m];
361: xs = 0;
362: for (i=0; i<(rank%m); i++) { xs += lx[i];}
363: if (m > 1 && x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %d %d",x,s);
364: } else if (flg2) {
365: SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
366: } else { /* Normal PETSc distribution */
367: x = M/m + ((M % m) > (rank % m));
368: if (m > 1 && x < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Column width is too thin for stencil! %d %d",x,s);
369: if ((M % m) > (rank % m)) { xs = (rank % m)*x; }
370: else { xs = (M % m)*(x+1) + ((rank % m)-(M % m))*x; }
371: PetscMalloc(m*sizeof(int),&lx);
372: flx = lx;
373: for (i=0; i<m; i++) {
374: lx[i] = M/m + ((M % m) > (i % m));
375: }
376: }
377: if (ly) { /* user decided distribution */
378: y = ly[(rank % (m*n))/m];
379: if (n > 1 && y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %d %d",y,s);
380: ys = 0;
381: for (i=0; i<(rank % (m*n))/m; i++) { ys += ly[i];}
382: } else if (flg2) {
383: SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
384: } else { /* Normal PETSc distribution */
385: y = N/n + ((N % n) > ((rank % (m*n)) /m));
386: if (n > 1 && y < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Row width is too thin for stencil! %d %d",y,s);
387: if ((N % n) > ((rank % (m*n)) /m)) {ys = ((rank % (m*n))/m)*y;}
388: else {ys = (N % n)*(y+1) + (((rank % (m*n))/m)-(N % n))*y;}
389: PetscMalloc(n*sizeof(int),&ly);
390: fly = ly;
391: for (i=0; i<n; i++) {
392: ly[i] = N/n + ((N % n) > (i % n));
393: }
394: }
395: if (lz) { /* user decided distribution */
396: z = lz[rank/(m*n)];
397: if (p > 1 && z < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Plane width is too thin for stencil! %d %d",z,s);
398: zs = 0;
399: for (i=0; i<(rank/(m*n)); i++) { zs += lz[i];}
400: } else if (flg2) {
401: SETERRQ(PETSC_ERR_SUP,"-da_partition_nodes_at_end not supported");
402: } else { /* Normal PETSc distribution */
403: z = P/p + ((P % p) > (rank / (m*n)));
404: if (p > 1 && z < s) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"Plane width is too thin for stencil! %d %d",z,s);
405: if ((P % p) > (rank / (m*n))) {zs = (rank/(m*n))*z;}
406: else {zs = (P % p)*(z+1) + ((rank/(m*n))-(P % p))*z;}
407: PetscMalloc(p*sizeof(int),&lz);
408: flz = lz;
409: for (i=0; i<p; i++) {
410: lz[i] = P/p + ((P % p) > (i % p));
411: }
412: }
413: ye = ys + y;
414: xe = xs + x;
415: ze = zs + z;
417: /* determine ghost region */
418: /* Assume No Periodicity */
419: if (xs-s > 0) Xs = xs - s; else Xs = 0;
420: if (ys-s > 0) Ys = ys - s; else Ys = 0;
421: if (zs-s > 0) Zs = zs - s; else Zs = 0;
422: if (xe+s <= M) Xe = xe + s; else Xe = M;
423: if (ye+s <= N) Ye = ye + s; else Ye = N;
424: if (ze+s <= P) Ze = ze + s; else Ze = P;
426: /* X Periodic */
427: if (DAXPeriodic(wrap)){
428: Xs = xs - s;
429: Xe = xe + s;
430: }
432: /* Y Periodic */
433: if (DAYPeriodic(wrap)){
434: Ys = ys - s;
435: Ye = ye + s;
436: }
438: /* Z Periodic */
439: if (DAZPeriodic(wrap)){
440: Zs = zs - s;
441: Ze = ze + s;
442: }
444: /* Resize all X parameters to reflect w */
445: x *= dof;
446: xs *= dof;
447: xe *= dof;
448: Xs *= dof;
449: Xe *= dof;
450: s_x = s*dof;
451: s_y = s;
452: s_z = s;
454: /* determine starting point of each processor */
455: nn = x*y*z;
456: PetscMalloc((2*size+1)*sizeof(int),&bases);
457: ldims = (int*)(bases+size+1);
458: MPI_Allgather(&nn,1,MPI_INT,ldims,1,MPI_INT,comm);
459: bases[0] = 0;
460: for (i=1; i<=size; i++) {
461: bases[i] = ldims[i-1];
462: }
463: for (i=1; i<=size; i++) {
464: bases[i] += bases[i-1];
465: }
467: /* allocate the base parallel and sequential vectors */
468: da->Nlocal = x*y*z;
469: VecCreateMPIWithArray(comm,da->Nlocal,PETSC_DECIDE,0,&global);
470: VecSetBlockSize(global,dof);
471: da->nlocal = (Xe-Xs)*(Ye-Ys)*(Ze-Zs);
472: VecCreateSeqWithArray(MPI_COMM_SELF,da->nlocal,0,&local);
473: VecSetBlockSize(local,dof);
475: /* generate appropriate vector scatters */
476: /* local to global inserts non-ghost point region into global */
477: VecGetOwnershipRange(global,&start,&end);
478: ISCreateStride(comm,x*y*z,start,1,&to);
480: left = xs - Xs;
481: bottom = ys - Ys; top = bottom + y;
482: down = zs - Zs; up = down + z;
483: count = x*(top-bottom)*(up-down);
484: PetscMalloc(count*sizeof(int),&idx);
485: count = 0;
486: for (i=down; i<up; i++) {
487: for (j=bottom; j<top; j++) {
488: for (k=0; k<x; k++) {
489: idx[count++] = (left+j*(Xe-Xs))+i*(Xe-Xs)*(Ye-Ys) + k;
490: }
491: }
492: }
493: ISCreateGeneral(comm,count,idx,&from);
494: PetscFree(idx);
496: VecScatterCreate(local,from,global,to,<og);
497: PetscLogObjectParent(da,to);
498: PetscLogObjectParent(da,from);
499: PetscLogObjectParent(da,ltog);
500: ISDestroy(from);
501: ISDestroy(to);
503: /* global to local must include ghost points */
504: if (stencil_type == DA_STENCIL_BOX) {
505: ISCreateStride(comm,(Xe-Xs)*(Ye-Ys)*(Ze-Zs),0,1,&to);
506: } else {
507: /* This is way ugly! We need to list the funny cross type region */
508: /* the bottom chunck */
509: left = xs - Xs;
510: bottom = ys - Ys; top = bottom + y;
511: down = zs - Zs; up = down + z;
512: count = down*(top-bottom)*x +
513: (up-down)*(bottom*x + (top-bottom)*(Xe-Xs) + (Ye-Ys-top)*x) +
514: (Ze-Zs-up)*(top-bottom)*x;
515: PetscMalloc(count*sizeof(int),&idx);
516: count = 0;
517: for (i=0; i<down; i++) {
518: for (j=bottom; j<top; j++) {
519: for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
520: }
521: }
522: /* the middle piece */
523: for (i=down; i<up; i++) {
524: /* front */
525: for (j=0; j<bottom; j++) {
526: for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
527: }
528: /* middle */
529: for (j=bottom; j<top; j++) {
530: for (k=0; k<Xe-Xs; k++) idx[count++] = j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
531: }
532: /* back */
533: for (j=top; j<Ye-Ys; j++) {
534: for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
535: }
536: }
537: /* the top piece */
538: for (i=up; i<Ze-Zs; i++) {
539: for (j=bottom; j<top; j++) {
540: for (k=0; k<x; k++) idx[count++] = left+j*(Xe-Xs)+i*(Xe-Xs)*(Ye-Ys)+k;
541: }
542: }
543: ISCreateGeneral(comm,count,idx,&to);
544: PetscFree(idx);
545: }
547: /* determine who lies on each side of use stored in n24 n25 n26
548: n21 n22 n23
549: n18 n19 n20
551: n15 n16 n17
552: n12 n14
553: n9 n10 n11
555: n6 n7 n8
556: n3 n4 n5
557: n0 n1 n2
558: */
559:
560: /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
561:
562: /* Assume Nodes are Internal to the Cube */
563:
564: n0 = rank - m*n - m - 1;
565: n1 = rank - m*n - m;
566: n2 = rank - m*n - m + 1;
567: n3 = rank - m*n -1;
568: n4 = rank - m*n;
569: n5 = rank - m*n + 1;
570: n6 = rank - m*n + m - 1;
571: n7 = rank - m*n + m;
572: n8 = rank - m*n + m + 1;
574: n9 = rank - m - 1;
575: n10 = rank - m;
576: n11 = rank - m + 1;
577: n12 = rank - 1;
578: n14 = rank + 1;
579: n15 = rank + m - 1;
580: n16 = rank + m;
581: n17 = rank + m + 1;
583: n18 = rank + m*n - m - 1;
584: n19 = rank + m*n - m;
585: n20 = rank + m*n - m + 1;
586: n21 = rank + m*n - 1;
587: n22 = rank + m*n;
588: n23 = rank + m*n + 1;
589: n24 = rank + m*n + m - 1;
590: n25 = rank + m*n + m;
591: n26 = rank + m*n + m + 1;
593: /* Assume Pieces are on Faces of Cube */
595: if (xs == 0) { /* First assume not corner or edge */
596: n0 = rank -1 - (m*n);
597: n3 = rank + m -1 - (m*n);
598: n6 = rank + 2*m -1 - (m*n);
599: n9 = rank -1;
600: n12 = rank + m -1;
601: n15 = rank + 2*m -1;
602: n18 = rank -1 + (m*n);
603: n21 = rank + m -1 + (m*n);
604: n24 = rank + 2*m -1 + (m*n);
605: }
607: if (xe == M*dof) { /* First assume not corner or edge */
608: n2 = rank -2*m +1 - (m*n);
609: n5 = rank - m +1 - (m*n);
610: n8 = rank +1 - (m*n);
611: n11 = rank -2*m +1;
612: n14 = rank - m +1;
613: n17 = rank +1;
614: n20 = rank -2*m +1 + (m*n);
615: n23 = rank - m +1 + (m*n);
616: n26 = rank +1 + (m*n);
617: }
619: if (ys==0) { /* First assume not corner or edge */
620: n0 = rank + m * (n-1) -1 - (m*n);
621: n1 = rank + m * (n-1) - (m*n);
622: n2 = rank + m * (n-1) +1 - (m*n);
623: n9 = rank + m * (n-1) -1;
624: n10 = rank + m * (n-1);
625: n11 = rank + m * (n-1) +1;
626: n18 = rank + m * (n-1) -1 + (m*n);
627: n19 = rank + m * (n-1) + (m*n);
628: n20 = rank + m * (n-1) +1 + (m*n);
629: }
631: if (ye == N) { /* First assume not corner or edge */
632: n6 = rank - m * (n-1) -1 - (m*n);
633: n7 = rank - m * (n-1) - (m*n);
634: n8 = rank - m * (n-1) +1 - (m*n);
635: n15 = rank - m * (n-1) -1;
636: n16 = rank - m * (n-1);
637: n17 = rank - m * (n-1) +1;
638: n24 = rank - m * (n-1) -1 + (m*n);
639: n25 = rank - m * (n-1) + (m*n);
640: n26 = rank - m * (n-1) +1 + (m*n);
641: }
642:
643: if (zs == 0) { /* First assume not corner or edge */
644: n0 = size - (m*n) + rank - m - 1;
645: n1 = size - (m*n) + rank - m;
646: n2 = size - (m*n) + rank - m + 1;
647: n3 = size - (m*n) + rank - 1;
648: n4 = size - (m*n) + rank;
649: n5 = size - (m*n) + rank + 1;
650: n6 = size - (m*n) + rank + m - 1;
651: n7 = size - (m*n) + rank + m ;
652: n8 = size - (m*n) + rank + m + 1;
653: }
655: if (ze == P) { /* First assume not corner or edge */
656: n18 = (m*n) - (size-rank) - m - 1;
657: n19 = (m*n) - (size-rank) - m;
658: n20 = (m*n) - (size-rank) - m + 1;
659: n21 = (m*n) - (size-rank) - 1;
660: n22 = (m*n) - (size-rank);
661: n23 = (m*n) - (size-rank) + 1;
662: n24 = (m*n) - (size-rank) + m - 1;
663: n25 = (m*n) - (size-rank) + m;
664: n26 = (m*n) - (size-rank) + m + 1;
665: }
667: if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
668: n0 = size - m*n + rank + m-1 - m;
669: n3 = size - m*n + rank + m-1;
670: n6 = size - m*n + rank + m-1 + m;
671: }
672:
673: if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
674: n18 = m*n - (size - rank) + m-1 - m;
675: n21 = m*n - (size - rank) + m-1;
676: n24 = m*n - (size - rank) + m-1 + m;
677: }
679: if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
680: n0 = rank + m*n -1 - m*n;
681: n9 = rank + m*n -1;
682: n18 = rank + m*n -1 + m*n;
683: }
685: if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
686: n6 = rank - m*(n-1) + m-1 - m*n;
687: n15 = rank - m*(n-1) + m-1;
688: n24 = rank - m*(n-1) + m-1 + m*n;
689: }
691: if ((xe==M*dof) && (zs==0)) { /* Assume an edge, not corner */
692: n2 = size - (m*n-rank) - (m-1) - m;
693: n5 = size - (m*n-rank) - (m-1);
694: n8 = size - (m*n-rank) - (m-1) + m;
695: }
697: if ((xe==M*dof) && (ze==P)) { /* Assume an edge, not corner */
698: n20 = m*n - (size - rank) - (m-1) - m;
699: n23 = m*n - (size - rank) - (m-1);
700: n26 = m*n - (size - rank) - (m-1) + m;
701: }
703: if ((xe==M*dof) && (ys==0)) { /* Assume an edge, not corner */
704: n2 = rank + m*(n-1) - (m-1) - m*n;
705: n11 = rank + m*(n-1) - (m-1);
706: n20 = rank + m*(n-1) - (m-1) + m*n;
707: }
709: if ((xe==M*dof) && (ye==N)) { /* Assume an edge, not corner */
710: n8 = rank - m*n +1 - m*n;
711: n17 = rank - m*n +1;
712: n26 = rank - m*n +1 + m*n;
713: }
715: if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
716: n0 = size - m + rank -1;
717: n1 = size - m + rank;
718: n2 = size - m + rank +1;
719: }
721: if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
722: n18 = m*n - (size - rank) + m*(n-1) -1;
723: n19 = m*n - (size - rank) + m*(n-1);
724: n20 = m*n - (size - rank) + m*(n-1) +1;
725: }
727: if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
728: n6 = size - (m*n-rank) - m * (n-1) -1;
729: n7 = size - (m*n-rank) - m * (n-1);
730: n8 = size - (m*n-rank) - m * (n-1) +1;
731: }
733: if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
734: n24 = rank - (size-m) -1;
735: n25 = rank - (size-m);
736: n26 = rank - (size-m) +1;
737: }
739: /* Check for Corners */
740: if ((xs==0) && (ys==0) && (zs==0)) { n0 = size -1;}
741: if ((xs==0) && (ys==0) && (ze==P)) { n18 = m*n-1;}
742: if ((xs==0) && (ye==N) && (zs==0)) { n6 = (size-1)-m*(n-1);}
743: if ((xs==0) && (ye==N) && (ze==P)) { n24 = m-1;}
744: if ((xe==M*dof) && (ys==0) && (zs==0)) { n2 = size-m;}
745: if ((xe==M*dof) && (ys==0) && (ze==P)) { n20 = m*n-m;}
746: if ((xe==M*dof) && (ye==N) && (zs==0)) { n8 = size-m*n;}
747: if ((xe==M*dof) && (ye==N) && (ze==P)) { n26 = 0;}
749: /* Check for when not X,Y, and Z Periodic */
751: /* If not X periodic */
752: if ((wrap != DA_XPERIODIC) && (wrap != DA_XYPERIODIC) &&
753: (wrap != DA_XZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
754: if (xs==0) {n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2;}
755: if (xe==M*dof) {n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
756: }
758: /* If not Y periodic */
759: if ((wrap != DA_YPERIODIC) && (wrap != DA_XYPERIODIC) &&
760: (wrap != DA_YZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
761: if (ys==0) {n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2;}
762: if (ye==N) {n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
763: }
765: /* If not Z periodic */
766: if ((wrap != DA_ZPERIODIC) && (wrap != DA_XZPERIODIC) &&
767: (wrap != DA_YZPERIODIC) && (wrap != DA_XYZPERIODIC)) {
768: if (zs==0) {n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2;}
769: if (ze==P) {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
770: }
772: /* If star stencil then delete the corner neighbors */
773: if (stencil_type == DA_STENCIL_STAR) {
774: /* save information about corner neighbors */
775: sn0 = n0; sn1 = n1; sn2 = n2; sn3 = n3; sn5 = n5; sn6 = n6; sn7 = n7;
776: sn8 = n8; sn9 = n9; sn11 = n11; sn15 = n15; sn17 = n17; sn18 = n18;
777: sn19 = n19; sn20 = n20; sn21 = n21; sn23 = n23; sn24 = n24; sn25 = n25;
778: sn26 = n26;
779: n0 = n1 = n2 = n3 = n5 = n6 = n7 = n8 = n9 = n11 =
780: n15 = n17 = n18 = n19 = n20 = n21 = n23 = n24 = n25 = n26 = -1;
781: }
784: PetscMalloc((Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(int),&idx);
785: PetscLogObjectMemory(da,(Xe-Xs)*(Ye-Ys)*(Ze-Zs)*sizeof(int));
787: nn = 0;
789: /* Bottom Level */
790: for (k=0; k<s_z; k++) {
791: for (i=1; i<=s_y; i++) {
792: if (n0 >= 0) { /* left below */
793: x_t = lx[n0 % m]*dof;
794: y_t = ly[(n0 % (m*n))/m];
795: z_t = lz[n0 / (m*n)];
796: s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x - (s_z-k-1)*x_t*y_t;
797: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
798: }
799: if (n1 >= 0) { /* directly below */
800: x_t = x;
801: y_t = ly[(n1 % (m*n))/m];
802: z_t = lz[n1 / (m*n)];
803: s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
804: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
805: }
806: if (n2 >= 0) { /* right below */
807: x_t = lx[n2 % m]*dof;
808: y_t = ly[(n2 % (m*n))/m];
809: z_t = lz[n2 / (m*n)];
810: s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
811: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
812: }
813: }
815: for (i=0; i<y; i++) {
816: if (n3 >= 0) { /* directly left */
817: x_t = lx[n3 % m]*dof;
818: y_t = y;
819: z_t = lz[n3 / (m*n)];
820: s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
821: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
822: }
824: if (n4 >= 0) { /* middle */
825: x_t = x;
826: y_t = y;
827: z_t = lz[n4 / (m*n)];
828: s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
829: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
830: }
832: if (n5 >= 0) { /* directly right */
833: x_t = lx[n5 % m]*dof;
834: y_t = y;
835: z_t = lz[n5 / (m*n)];
836: s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
837: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
838: }
839: }
841: for (i=1; i<=s_y; i++) {
842: if (n6 >= 0) { /* left above */
843: x_t = lx[n6 % m]*dof;
844: y_t = ly[(n6 % (m*n))/m];
845: z_t = lz[n6 / (m*n)];
846: s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
847: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
848: }
849: if (n7 >= 0) { /* directly above */
850: x_t = x;
851: y_t = ly[(n7 % (m*n))/m];
852: z_t = lz[n7 / (m*n)];
853: s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
854: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
855: }
856: if (n8 >= 0) { /* right above */
857: x_t = lx[n8 % m]*dof;
858: y_t = ly[(n8 % (m*n))/m];
859: z_t = lz[n8 / (m*n)];
860: s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
861: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
862: }
863: }
864: }
866: /* Middle Level */
867: for (k=0; k<z; k++) {
868: for (i=1; i<=s_y; i++) {
869: if (n9 >= 0) { /* left below */
870: x_t = lx[n9 % m]*dof;
871: y_t = ly[(n9 % (m*n))/m];
872: /* z_t = z; */
873: s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
874: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
875: }
876: if (n10 >= 0) { /* directly below */
877: x_t = x;
878: y_t = ly[(n10 % (m*n))/m];
879: /* z_t = z; */
880: s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
881: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
882: }
883: if (n11 >= 0) { /* right below */
884: x_t = lx[n11 % m]*dof;
885: y_t = ly[(n11 % (m*n))/m];
886: /* z_t = z; */
887: s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
888: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
889: }
890: }
892: for (i=0; i<y; i++) {
893: if (n12 >= 0) { /* directly left */
894: x_t = lx[n12 % m]*dof;
895: y_t = y;
896: /* z_t = z; */
897: s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
898: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
899: }
901: /* Interior */
902: s_t = bases[rank] + i*x + k*x*y;
903: for (j=0; j<x; j++) { idx[nn++] = s_t++;}
905: if (n14 >= 0) { /* directly right */
906: x_t = lx[n14 % m]*dof;
907: y_t = y;
908: /* z_t = z; */
909: s_t = bases[n14] + i*x_t + k*x_t*y_t;
910: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
911: }
912: }
914: for (i=1; i<=s_y; i++) {
915: if (n15 >= 0) { /* left above */
916: x_t = lx[n15 % m]*dof;
917: y_t = ly[(n15 % (m*n))/m];
918: /* z_t = z; */
919: s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
920: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
921: }
922: if (n16 >= 0) { /* directly above */
923: x_t = x;
924: y_t = ly[(n16 % (m*n))/m];
925: /* z_t = z; */
926: s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
927: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
928: }
929: if (n17 >= 0) { /* right above */
930: x_t = lx[n17 % m]*dof;
931: y_t = ly[(n17 % (m*n))/m];
932: /* z_t = z; */
933: s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
934: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
935: }
936: }
937: }
938:
939: /* Upper Level */
940: for (k=0; k<s_z; k++) {
941: for (i=1; i<=s_y; i++) {
942: if (n18 >= 0) { /* left below */
943: x_t = lx[n18 % m]*dof;
944: y_t = ly[(n18 % (m*n))/m];
945: /* z_t = lz[n18 / (m*n)]; */
946: s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
947: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
948: }
949: if (n19 >= 0) { /* directly below */
950: x_t = x;
951: y_t = ly[(n19 % (m*n))/m];
952: /* z_t = lz[n19 / (m*n)]; */
953: s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
954: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
955: }
956: if (n20 >= 0) { /* right below */
957: x_t = lx[n20 % m]*dof;
958: y_t = ly[(n20 % (m*n))/m];
959: /* z_t = lz[n20 / (m*n)]; */
960: s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
961: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
962: }
963: }
965: for (i=0; i<y; i++) {
966: if (n21 >= 0) { /* directly left */
967: x_t = lx[n21 % m]*dof;
968: y_t = y;
969: /* z_t = lz[n21 / (m*n)]; */
970: s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
971: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
972: }
974: if (n22 >= 0) { /* middle */
975: x_t = x;
976: y_t = y;
977: /* z_t = lz[n22 / (m*n)]; */
978: s_t = bases[n22] + i*x_t + k*x_t*y_t;
979: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
980: }
982: if (n23 >= 0) { /* directly right */
983: x_t = lx[n23 % m]*dof;
984: y_t = y;
985: /* z_t = lz[n23 / (m*n)]; */
986: s_t = bases[n23] + i*x_t + k*x_t*y_t;
987: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
988: }
989: }
991: for (i=1; i<=s_y; i++) {
992: if (n24 >= 0) { /* left above */
993: x_t = lx[n24 % m]*dof;
994: y_t = ly[(n24 % (m*n))/m];
995: /* z_t = lz[n24 / (m*n)]; */
996: s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
997: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
998: }
999: if (n25 >= 0) { /* directly above */
1000: x_t = x;
1001: y_t = ly[(n25 % (m*n))/m];
1002: /* z_t = lz[n25 / (m*n)]; */
1003: s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1004: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1005: }
1006: if (n26 >= 0) { /* right above */
1007: x_t = lx[n26 % m]*dof;
1008: y_t = ly[(n26 % (m*n))/m];
1009: /* z_t = lz[n26 / (m*n)]; */
1010: s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1011: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1012: }
1013: }
1014: }
1015: base = bases[rank];
1016: ISCreateGeneral(comm,nn,idx,&from);
1017: VecScatterCreate(global,from,local,to,>ol);
1018: PetscLogObjectParent(da,gtol);
1019: PetscLogObjectParent(da,to);
1020: PetscLogObjectParent(da,from);
1021: ISDestroy(to);
1022: ISDestroy(from);
1023: da->stencil_type = stencil_type;
1024: da->M = M; da->N = N; da->P = P;
1025: da->m = m; da->n = n; da->p = p;
1026: da->w = dof; da->s = s;
1027: da->xs = xs; da->xe = xe; da->ys = ys; da->ye = ye; da->zs = zs; da->ze = ze;
1028: da->Xs = Xs; da->Xe = Xe; da->Ys = Ys; da->Ye = Ye; da->Zs = Zs; da->Ze = Ze;
1030: VecDestroy(local);
1031: VecDestroy(global);
1033: if (stencil_type == DA_STENCIL_STAR) {
1034: /*
1035: Recompute the local to global mappings, this time keeping the
1036: information about the cross corner processor numbers.
1037: */
1038: n0 = sn0; n1 = sn1; n2 = sn2; n3 = sn3; n5 = sn5; n6 = sn6; n7 = sn7;
1039: n8 = sn8; n9 = sn9; n11 = sn11; n15 = sn15; n17 = sn17; n18 = sn18;
1040: n19 = sn19; n20 = sn20; n21 = sn21; n23 = sn23; n24 = sn24; n25 = sn25;
1041: n26 = sn26;
1043: nn = 0;
1045: /* Bottom Level */
1046: for (k=0; k<s_z; k++) {
1047: for (i=1; i<=s_y; i++) {
1048: if (n0 >= 0) { /* left below */
1049: x_t = lx[n0 % m]*dof;
1050: y_t = ly[(n0 % (m*n))/m];
1051: z_t = lz[n0 / (m*n)];
1052: s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t - s_x - (s_z-k-1)*x_t*y_t;
1053: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1054: }
1055: if (n1 >= 0) { /* directly below */
1056: x_t = x;
1057: y_t = ly[(n1 % (m*n))/m];
1058: z_t = lz[n1 / (m*n)];
1059: s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1060: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1061: }
1062: if (n2 >= 0) { /* right below */
1063: x_t = lx[n2 % m]*dof;
1064: y_t = ly[(n2 % (m*n))/m];
1065: z_t = lz[n2 / (m*n)];
1066: s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1067: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1068: }
1069: }
1071: for (i=0; i<y; i++) {
1072: if (n3 >= 0) { /* directly left */
1073: x_t = lx[n3 % m]*dof;
1074: y_t = y;
1075: z_t = lz[n3 / (m*n)];
1076: s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1077: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1078: }
1080: if (n4 >= 0) { /* middle */
1081: x_t = x;
1082: y_t = y;
1083: z_t = lz[n4 / (m*n)];
1084: s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1085: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1086: }
1088: if (n5 >= 0) { /* directly right */
1089: x_t = lx[n5 % m]*dof;
1090: y_t = y;
1091: z_t = lz[n5 / (m*n)];
1092: s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1093: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1094: }
1095: }
1097: for (i=1; i<=s_y; i++) {
1098: if (n6 >= 0) { /* left above */
1099: x_t = lx[n6 % m]*dof;
1100: y_t = ly[(n6 % (m*n))/m];
1101: z_t = lz[n6 / (m*n)];
1102: s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1103: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1104: }
1105: if (n7 >= 0) { /* directly above */
1106: x_t = x;
1107: y_t = ly[(n7 % (m*n))/m];
1108: z_t = lz[n7 / (m*n)];
1109: s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1110: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1111: }
1112: if (n8 >= 0) { /* right above */
1113: x_t = lx[n8 % m]*dof;
1114: y_t = ly[(n8 % (m*n))/m];
1115: z_t = lz[n8 / (m*n)];
1116: s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1117: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1118: }
1119: }
1120: }
1122: /* Middle Level */
1123: for (k=0; k<z; k++) {
1124: for (i=1; i<=s_y; i++) {
1125: if (n9 >= 0) { /* left below */
1126: x_t = lx[n9 % m]*dof;
1127: y_t = ly[(n9 % (m*n))/m];
1128: /* z_t = z; */
1129: s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1130: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1131: }
1132: if (n10 >= 0) { /* directly below */
1133: x_t = x;
1134: y_t = ly[(n10 % (m*n))/m];
1135: /* z_t = z; */
1136: s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1137: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1138: }
1139: if (n11 >= 0) { /* right below */
1140: x_t = lx[n11 % m]*dof;
1141: y_t = ly[(n11 % (m*n))/m];
1142: /* z_t = z; */
1143: s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1144: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1145: }
1146: }
1148: for (i=0; i<y; i++) {
1149: if (n12 >= 0) { /* directly left */
1150: x_t = lx[n12 % m]*dof;
1151: y_t = y;
1152: /* z_t = z; */
1153: s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
1154: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1155: }
1157: /* Interior */
1158: s_t = bases[rank] + i*x + k*x*y;
1159: for (j=0; j<x; j++) { idx[nn++] = s_t++;}
1161: if (n14 >= 0) { /* directly right */
1162: x_t = lx[n14 % m]*dof;
1163: y_t = y;
1164: /* z_t = z; */
1165: s_t = bases[n14] + i*x_t + k*x_t*y_t;
1166: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1167: }
1168: }
1170: for (i=1; i<=s_y; i++) {
1171: if (n15 >= 0) { /* left above */
1172: x_t = lx[n15 % m]*dof;
1173: y_t = ly[(n15 % (m*n))/m];
1174: /* z_t = z; */
1175: s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
1176: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1177: }
1178: if (n16 >= 0) { /* directly above */
1179: x_t = x;
1180: y_t = ly[(n16 % (m*n))/m];
1181: /* z_t = z; */
1182: s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
1183: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1184: }
1185: if (n17 >= 0) { /* right above */
1186: x_t = lx[n17 % m]*dof;
1187: y_t = ly[(n17 % (m*n))/m];
1188: /* z_t = z; */
1189: s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
1190: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1191: }
1192: }
1193: }
1194:
1195: /* Upper Level */
1196: for (k=0; k<s_z; k++) {
1197: for (i=1; i<=s_y; i++) {
1198: if (n18 >= 0) { /* left below */
1199: x_t = lx[n18 % m]*dof;
1200: y_t = ly[(n18 % (m*n))/m];
1201: /* z_t = lz[n18 / (m*n)]; */
1202: s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1203: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1204: }
1205: if (n19 >= 0) { /* directly below */
1206: x_t = x;
1207: y_t = ly[(n19 % (m*n))/m];
1208: /* z_t = lz[n19 / (m*n)]; */
1209: s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1210: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1211: }
1212: if (n20 >= 0) { /* right below */
1213: x_t = lx[n20 % m]*dof;
1214: y_t = ly[(n20 % (m*n))/m];
1215: /* z_t = lz[n20 / (m*n)]; */
1216: s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1217: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1218: }
1219: }
1221: for (i=0; i<y; i++) {
1222: if (n21 >= 0) { /* directly left */
1223: x_t = lx[n21 % m]*dof;
1224: y_t = y;
1225: /* z_t = lz[n21 / (m*n)]; */
1226: s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
1227: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1228: }
1230: if (n22 >= 0) { /* middle */
1231: x_t = x;
1232: y_t = y;
1233: /* z_t = lz[n22 / (m*n)]; */
1234: s_t = bases[n22] + i*x_t + k*x_t*y_t;
1235: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1236: }
1238: if (n23 >= 0) { /* directly right */
1239: x_t = lx[n23 % m]*dof;
1240: y_t = y;
1241: /* z_t = lz[n23 / (m*n)]; */
1242: s_t = bases[n23] + i*x_t + k*x_t*y_t;
1243: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1244: }
1245: }
1247: for (i=1; i<=s_y; i++) {
1248: if (n24 >= 0) { /* left above */
1249: x_t = lx[n24 % m]*dof;
1250: y_t = ly[(n24 % (m*n))/m];
1251: /* z_t = lz[n24 / (m*n)]; */
1252: s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
1253: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1254: }
1255: if (n25 >= 0) { /* directly above */
1256: x_t = x;
1257: y_t = ly[(n25 % (m*n))/m];
1258: /* z_t = lz[n25 / (m*n)]; */
1259: s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1260: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1261: }
1262: if (n26 >= 0) { /* right above */
1263: x_t = lx[n26 % m]*dof;
1264: y_t = ly[(n26 % (m*n))/m];
1265: /* z_t = lz[n26 / (m*n)]; */
1266: s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1267: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1268: }
1269: }
1270: }
1271: }
1272: /* redo idx to include "missing" ghost points */
1273: /* Solve for X,Y, and Z Periodic Case First, Then Modify Solution */
1274:
1275: /* Assume Nodes are Internal to the Cube */
1276:
1277: n0 = rank - m*n - m - 1;
1278: n1 = rank - m*n - m;
1279: n2 = rank - m*n - m + 1;
1280: n3 = rank - m*n -1;
1281: n4 = rank - m*n;
1282: n5 = rank - m*n + 1;
1283: n6 = rank - m*n + m - 1;
1284: n7 = rank - m*n + m;
1285: n8 = rank - m*n + m + 1;
1287: n9 = rank - m - 1;
1288: n10 = rank - m;
1289: n11 = rank - m + 1;
1290: n12 = rank - 1;
1291: n14 = rank + 1;
1292: n15 = rank + m - 1;
1293: n16 = rank + m;
1294: n17 = rank + m + 1;
1296: n18 = rank + m*n - m - 1;
1297: n19 = rank + m*n - m;
1298: n20 = rank + m*n - m + 1;
1299: n21 = rank + m*n - 1;
1300: n22 = rank + m*n;
1301: n23 = rank + m*n + 1;
1302: n24 = rank + m*n + m - 1;
1303: n25 = rank + m*n + m;
1304: n26 = rank + m*n + m + 1;
1306: /* Assume Pieces are on Faces of Cube */
1308: if (xs == 0) { /* First assume not corner or edge */
1309: n0 = rank -1 - (m*n);
1310: n3 = rank + m -1 - (m*n);
1311: n6 = rank + 2*m -1 - (m*n);
1312: n9 = rank -1;
1313: n12 = rank + m -1;
1314: n15 = rank + 2*m -1;
1315: n18 = rank -1 + (m*n);
1316: n21 = rank + m -1 + (m*n);
1317: n24 = rank + 2*m -1 + (m*n);
1318: }
1320: if (xe == M*dof) { /* First assume not corner or edge */
1321: n2 = rank -2*m +1 - (m*n);
1322: n5 = rank - m +1 - (m*n);
1323: n8 = rank +1 - (m*n);
1324: n11 = rank -2*m +1;
1325: n14 = rank - m +1;
1326: n17 = rank +1;
1327: n20 = rank -2*m +1 + (m*n);
1328: n23 = rank - m +1 + (m*n);
1329: n26 = rank +1 + (m*n);
1330: }
1332: if (ys==0) { /* First assume not corner or edge */
1333: n0 = rank + m * (n-1) -1 - (m*n);
1334: n1 = rank + m * (n-1) - (m*n);
1335: n2 = rank + m * (n-1) +1 - (m*n);
1336: n9 = rank + m * (n-1) -1;
1337: n10 = rank + m * (n-1);
1338: n11 = rank + m * (n-1) +1;
1339: n18 = rank + m * (n-1) -1 + (m*n);
1340: n19 = rank + m * (n-1) + (m*n);
1341: n20 = rank + m * (n-1) +1 + (m*n);
1342: }
1344: if (ye == N) { /* First assume not corner or edge */
1345: n6 = rank - m * (n-1) -1 - (m*n);
1346: n7 = rank - m * (n-1) - (m*n);
1347: n8 = rank - m * (n-1) +1 - (m*n);
1348: n15 = rank - m * (n-1) -1;
1349: n16 = rank - m * (n-1);
1350: n17 = rank - m * (n-1) +1;
1351: n24 = rank - m * (n-1) -1 + (m*n);
1352: n25 = rank - m * (n-1) + (m*n);
1353: n26 = rank - m * (n-1) +1 + (m*n);
1354: }
1355:
1356: if (zs == 0) { /* First assume not corner or edge */
1357: n0 = size - (m*n) + rank - m - 1;
1358: n1 = size - (m*n) + rank - m;
1359: n2 = size - (m*n) + rank - m + 1;
1360: n3 = size - (m*n) + rank - 1;
1361: n4 = size - (m*n) + rank;
1362: n5 = size - (m*n) + rank + 1;
1363: n6 = size - (m*n) + rank + m - 1;
1364: n7 = size - (m*n) + rank + m ;
1365: n8 = size - (m*n) + rank + m + 1;
1366: }
1368: if (ze == P) { /* First assume not corner or edge */
1369: n18 = (m*n) - (size-rank) - m - 1;
1370: n19 = (m*n) - (size-rank) - m;
1371: n20 = (m*n) - (size-rank) - m + 1;
1372: n21 = (m*n) - (size-rank) - 1;
1373: n22 = (m*n) - (size-rank);
1374: n23 = (m*n) - (size-rank) + 1;
1375: n24 = (m*n) - (size-rank) + m - 1;
1376: n25 = (m*n) - (size-rank) + m;
1377: n26 = (m*n) - (size-rank) + m + 1;
1378: }
1380: if ((xs==0) && (zs==0)) { /* Assume an edge, not corner */
1381: n0 = size - m*n + rank + m-1 - m;
1382: n3 = size - m*n + rank + m-1;
1383: n6 = size - m*n + rank + m-1 + m;
1384: }
1385:
1386: if ((xs==0) && (ze==P)) { /* Assume an edge, not corner */
1387: n18 = m*n - (size - rank) + m-1 - m;
1388: n21 = m*n - (size - rank) + m-1;
1389: n24 = m*n - (size - rank) + m-1 + m;
1390: }
1392: if ((xs==0) && (ys==0)) { /* Assume an edge, not corner */
1393: n0 = rank + m*n -1 - m*n;
1394: n9 = rank + m*n -1;
1395: n18 = rank + m*n -1 + m*n;
1396: }
1398: if ((xs==0) && (ye==N)) { /* Assume an edge, not corner */
1399: n6 = rank - m*(n-1) + m-1 - m*n;
1400: n15 = rank - m*(n-1) + m-1;
1401: n24 = rank - m*(n-1) + m-1 + m*n;
1402: }
1404: if ((xe==M*dof) && (zs==0)) { /* Assume an edge, not corner */
1405: n2 = size - (m*n-rank) - (m-1) - m;
1406: n5 = size - (m*n-rank) - (m-1);
1407: n8 = size - (m*n-rank) - (m-1) + m;
1408: }
1410: if ((xe==M*dof) && (ze==P)) { /* Assume an edge, not corner */
1411: n20 = m*n - (size - rank) - (m-1) - m;
1412: n23 = m*n - (size - rank) - (m-1);
1413: n26 = m*n - (size - rank) - (m-1) + m;
1414: }
1416: if ((xe==M*dof) && (ys==0)) { /* Assume an edge, not corner */
1417: n2 = rank + m*(n-1) - (m-1) - m*n;
1418: n11 = rank + m*(n-1) - (m-1);
1419: n20 = rank + m*(n-1) - (m-1) + m*n;
1420: }
1422: if ((xe==M*dof) && (ye==N)) { /* Assume an edge, not corner */
1423: n8 = rank - m*n +1 - m*n;
1424: n17 = rank - m*n +1;
1425: n26 = rank - m*n +1 + m*n;
1426: }
1428: if ((ys==0) && (zs==0)) { /* Assume an edge, not corner */
1429: n0 = size - m + rank -1;
1430: n1 = size - m + rank;
1431: n2 = size - m + rank +1;
1432: }
1434: if ((ys==0) && (ze==P)) { /* Assume an edge, not corner */
1435: n18 = m*n - (size - rank) + m*(n-1) -1;
1436: n19 = m*n - (size - rank) + m*(n-1);
1437: n20 = m*n - (size - rank) + m*(n-1) +1;
1438: }
1440: if ((ye==N) && (zs==0)) { /* Assume an edge, not corner */
1441: n6 = size - (m*n-rank) - m * (n-1) -1;
1442: n7 = size - (m*n-rank) - m * (n-1);
1443: n8 = size - (m*n-rank) - m * (n-1) +1;
1444: }
1446: if ((ye==N) && (ze==P)) { /* Assume an edge, not corner */
1447: n24 = rank - (size-m) -1;
1448: n25 = rank - (size-m);
1449: n26 = rank - (size-m) +1;
1450: }
1452: /* Check for Corners */
1453: if ((xs==0) && (ys==0) && (zs==0)) { n0 = size -1;}
1454: if ((xs==0) && (ys==0) && (ze==P)) { n18 = m*n-1;}
1455: if ((xs==0) && (ye==N) && (zs==0)) { n6 = (size-1)-m*(n-1);}
1456: if ((xs==0) && (ye==N) && (ze==P)) { n24 = m-1;}
1457: if ((xe==M*dof) && (ys==0) && (zs==0)) { n2 = size-m;}
1458: if ((xe==M*dof) && (ys==0) && (ze==P)) { n20 = m*n-m;}
1459: if ((xe==M*dof) && (ye==N) && (zs==0)) { n8 = size-m*n;}
1460: if ((xe==M*dof) && (ye==N) && (ze==P)) { n26 = 0;}
1462: /* Check for when not X,Y, and Z Periodic */
1464: /* If not X periodic */
1465: if (!DAXPeriodic(wrap)){
1466: if (xs==0) {n0 = n3 = n6 = n9 = n12 = n15 = n18 = n21 = n24 = -2;}
1467: if (xe==M*dof) {n2 = n5 = n8 = n11 = n14 = n17 = n20 = n23 = n26 = -2;}
1468: }
1470: /* If not Y periodic */
1471: if (!DAYPeriodic(wrap)){
1472: if (ys==0) {n0 = n1 = n2 = n9 = n10 = n11 = n18 = n19 = n20 = -2;}
1473: if (ye==N) {n6 = n7 = n8 = n15 = n16 = n17 = n24 = n25 = n26 = -2;}
1474: }
1476: /* If not Z periodic */
1477: if (!DAZPeriodic(wrap)){
1478: if (zs==0) {n0 = n1 = n2 = n3 = n4 = n5 = n6 = n7 = n8 = -2;}
1479: if (ze==P) {n18 = n19 = n20 = n21 = n22 = n23 = n24 = n25 = n26 = -2;}
1480: }
1482: nn = 0;
1484: /* Bottom Level */
1485: for (k=0; k<s_z; k++) {
1486: for (i=1; i<=s_y; i++) {
1487: if (n0 >= 0) { /* left below */
1488: x_t = lx[n0 % m]*dof;
1489: y_t = ly[(n0 % (m*n))/m];
1490: z_t = lz[n0 / (m*n)];
1491: s_t = bases[n0] + x_t*y_t*z_t - (s_y-i)*x_t -s_x - (s_z-k-1)*x_t*y_t;
1492: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1493: }
1494: if (n1 >= 0) { /* directly below */
1495: x_t = x;
1496: y_t = ly[(n1 % (m*n))/m];
1497: z_t = lz[n1 / (m*n)];
1498: s_t = bases[n1] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1499: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1500: }
1501: if (n2 >= 0) { /* right below */
1502: x_t = lx[n2 % m]*dof;
1503: y_t = ly[(n2 % (m*n))/m];
1504: z_t = lz[n2 / (m*n)];
1505: s_t = bases[n2] + x_t*y_t*z_t - (s_y+1-i)*x_t - (s_z-k-1)*x_t*y_t;
1506: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1507: }
1508: }
1510: for (i=0; i<y; i++) {
1511: if (n3 >= 0) { /* directly left */
1512: x_t = lx[n3 % m]*dof;
1513: y_t = y;
1514: z_t = lz[n3 / (m*n)];
1515: s_t = bases[n3] + (i+1)*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1516: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1517: }
1519: if (n4 >= 0) { /* middle */
1520: x_t = x;
1521: y_t = y;
1522: z_t = lz[n4 / (m*n)];
1523: s_t = bases[n4] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1524: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1525: }
1527: if (n5 >= 0) { /* directly right */
1528: x_t = lx[n5 % m]*dof;
1529: y_t = y;
1530: z_t = lz[n5 / (m*n)];
1531: s_t = bases[n5] + i*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1532: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1533: }
1534: }
1536: for (i=1; i<=s_y; i++) {
1537: if (n6 >= 0) { /* left above */
1538: x_t = lx[n6 % m]*dof;
1539: y_t = ly[(n6 % (m*n))/m];
1540: z_t = lz[n6 / (m*n)];
1541: s_t = bases[n6] + i*x_t - s_x + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1542: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1543: }
1544: if (n7 >= 0) { /* directly above */
1545: x_t = x;
1546: y_t = ly[(n7 % (m*n))/m];
1547: z_t = lz[n7 / (m*n)];
1548: s_t = bases[n7] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1549: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1550: }
1551: if (n8 >= 0) { /* right above */
1552: x_t = lx[n8 % m]*dof;
1553: y_t = ly[(n8 % (m*n))/m];
1554: z_t = lz[n8 / (m*n)];
1555: s_t = bases[n8] + (i-1)*x_t + x_t*y_t*z_t - (s_z-k)*x_t*y_t;
1556: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1557: }
1558: }
1559: }
1561: /* Middle Level */
1562: for (k=0; k<z; k++) {
1563: for (i=1; i<=s_y; i++) {
1564: if (n9 >= 0) { /* left below */
1565: x_t = lx[n9 % m]*dof;
1566: y_t = ly[(n9 % (m*n))/m];
1567: /* z_t = z; */
1568: s_t = bases[n9] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1569: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1570: }
1571: if (n10 >= 0) { /* directly below */
1572: x_t = x;
1573: y_t = ly[(n10 % (m*n))/m];
1574: /* z_t = z; */
1575: s_t = bases[n10] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1576: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1577: }
1578: if (n11 >= 0) { /* right below */
1579: x_t = lx[n11 % m]*dof;
1580: y_t = ly[(n11 % (m*n))/m];
1581: /* z_t = z; */
1582: s_t = bases[n11] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1583: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1584: }
1585: }
1587: for (i=0; i<y; i++) {
1588: if (n12 >= 0) { /* directly left */
1589: x_t = lx[n12 % m]*dof;
1590: y_t = y;
1591: /* z_t = z; */
1592: s_t = bases[n12] + (i+1)*x_t - s_x + k*x_t*y_t;
1593: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1594: }
1596: /* Interior */
1597: s_t = bases[rank] + i*x + k*x*y;
1598: for (j=0; j<x; j++) { idx[nn++] = s_t++;}
1600: if (n14 >= 0) { /* directly right */
1601: x_t = lx[n14 % m]*dof;
1602: y_t = y;
1603: /* z_t = z; */
1604: s_t = bases[n14] + i*x_t + k*x_t*y_t;
1605: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1606: }
1607: }
1609: for (i=1; i<=s_y; i++) {
1610: if (n15 >= 0) { /* left above */
1611: x_t = lx[n15 % m]*dof;
1612: y_t = ly[(n15 % (m*n))/m];
1613: /* z_t = z; */
1614: s_t = bases[n15] + i*x_t - s_x + k*x_t*y_t;
1615: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1616: }
1617: if (n16 >= 0) { /* directly above */
1618: x_t = x;
1619: y_t = ly[(n16 % (m*n))/m];
1620: /* z_t = z; */
1621: s_t = bases[n16] + (i-1)*x_t + k*x_t*y_t;
1622: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1623: }
1624: if (n17 >= 0) { /* right above */
1625: x_t = lx[n17 % m]*dof;
1626: y_t = ly[(n17 % (m*n))/m];
1627: /* z_t = z; */
1628: s_t = bases[n17] + (i-1)*x_t + k*x_t*y_t;
1629: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1630: }
1631: }
1632: }
1633:
1634: /* Upper Level */
1635: for (k=0; k<s_z; k++) {
1636: for (i=1; i<=s_y; i++) {
1637: if (n18 >= 0) { /* left below */
1638: x_t = lx[n18 % m]*dof;
1639: y_t = ly[(n18 % (m*n))/m];
1640: /* z_t = lz[n18 / (m*n)]; */
1641: s_t = bases[n18] - (s_y-i)*x_t -s_x + (k+1)*x_t*y_t;
1642: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1643: }
1644: if (n19 >= 0) { /* directly below */
1645: x_t = x;
1646: y_t = ly[(n19 % (m*n))/m];
1647: /* z_t = lz[n19 / (m*n)]; */
1648: s_t = bases[n19] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1649: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1650: }
1651: if (n20 >= 0) { /* right belodof */
1652: x_t = lx[n20 % m]*dof;
1653: y_t = ly[(n20 % (m*n))/m];
1654: /* z_t = lz[n20 / (m*n)]; */
1655: s_t = bases[n20] - (s_y+1-i)*x_t + (k+1)*x_t*y_t;
1656: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1657: }
1658: }
1660: for (i=0; i<y; i++) {
1661: if (n21 >= 0) { /* directly left */
1662: x_t = lx[n21 % m]*dof;
1663: y_t = y;
1664: /* z_t = lz[n21 / (m*n)]; */
1665: s_t = bases[n21] + (i+1)*x_t - s_x + k*x_t*y_t;
1666: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1667: }
1669: if (n22 >= 0) { /* middle */
1670: x_t = x;
1671: y_t = y;
1672: /* z_t = lz[n22 / (m*n)]; */
1673: s_t = bases[n22] + i*x_t + k*x_t*y_t;
1674: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1675: }
1677: if (n23 >= 0) { /* directly right */
1678: x_t = lx[n23 % m]*dof;
1679: y_t = y;
1680: /* z_t = lz[n23 / (m*n)]; */
1681: s_t = bases[n23] + i*x_t + k*x_t*y_t;
1682: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1683: }
1684: }
1686: for (i=1; i<=s_y; i++) {
1687: if (n24 >= 0) { /* left above */
1688: x_t = lx[n24 % m]*dof;
1689: y_t = ly[(n24 % (m*n))/m];
1690: /* z_t = lz[n24 / (m*n)]; */
1691: s_t = bases[n24] + i*x_t - s_x + k*x_t*y_t;
1692: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1693: }
1694: if (n25 >= 0) { /* directly above */
1695: x_t = x;
1696: y_t = ly[(n25 % (m*n))/m];
1697: /* z_t = lz[n25 / (m*n)]; */
1698: s_t = bases[n25] + (i-1)*x_t + k*x_t*y_t;
1699: for (j=0; j<x_t; j++) { idx[nn++] = s_t++;}
1700: }
1701: if (n26 >= 0) { /* right above */
1702: x_t = lx[n26 % m]*dof;
1703: y_t = ly[(n26 % (m*n))/m];
1704: /* z_t = lz[n26 / (m*n)]; */
1705: s_t = bases[n26] + (i-1)*x_t + k*x_t*y_t;
1706: for (j=0; j<s_x; j++) { idx[nn++] = s_t++;}
1707: }
1708: }
1709: }
1710: PetscFree(bases);
1711: da->gtol = gtol;
1712: da->ltog = ltog;
1713: da->idx = idx;
1714: da->Nl = nn;
1715: da->base = base;
1716: da->ops->view = DAView_3d;
1717: da->wrap = wrap;
1718: *inra = da;
1720: /*
1721: Set the local to global ordering in the global vector, this allows use
1722: of VecSetValuesLocal().
1723: */
1724: ISLocalToGlobalMappingCreateNC(comm,nn,idx,&da->ltogmap);
1725: ISLocalToGlobalMappingBlock(da->ltogmap,da->w,&da->ltogmapb);
1726: PetscLogObjectParent(da,da->ltogmap);
1728: da->ltol = PETSC_NULL;
1729: da->ao = PETSC_NULL;
1731: if (!flx) {
1732: PetscMalloc(m*sizeof(int),&flx);
1733: PetscMemcpy(flx,lx,m*sizeof(int));
1734: }
1735: if (!fly) {
1736: PetscMalloc(n*sizeof(int),&fly);
1737: PetscMemcpy(fly,ly,n*sizeof(int));
1738: }
1739: if (!flz) {
1740: PetscMalloc(p*sizeof(int),&flz);
1741: PetscMemcpy(flz,lz,p*sizeof(int));
1742: }
1743: da->lx = flx;
1744: da->ly = fly;
1745: da->lz = flz;
1747: PetscOptionsHasName(PETSC_NULL,"-da_view",&flg1);
1748: if (flg1) {DAView(da,PETSC_VIEWER_STDOUT_(da->comm));}
1749: PetscOptionsHasName(PETSC_NULL,"-da_view_draw",&flg1);
1750: if (flg1) {DAView(da,PETSC_VIEWER_DRAW_(da->comm));}
1751: PetscOptionsHasName(PETSC_NULL,"-help",&flg1);
1752: if (flg1) {DAPrintHelp(da);}
1753: PetscPublishAll(da);
1755: return(0);
1756: }