Actual source code: stagda.c
petsc-3.14.0 2020-09-29
1: /* Routines to convert between a (subset of) DMStag and DMDA */
3: #include <petscdmda.h>
4: #include <petsc/private/dmstagimpl.h>
5: #include <petscdmdatypes.h>
7: static PetscErrorCode DMStagCreateCompatibleDMDA(DM dm,DMStagStencilLocation loc,PetscInt c,DM *dmda)
8: {
9: PetscErrorCode ierr;
10: DM_Stag * const stag = (DM_Stag*) dm->data;
11: PetscInt dim,i,j,stencilWidth,dof,N[DMSTAG_MAX_DIM];
12: DMDAStencilType stencilType;
13: PetscInt *l[DMSTAG_MAX_DIM];
17: DMGetDimension(dm,&dim);
19: /* Create grid decomposition (to be adjusted later) */
20: for (i=0; i<dim; ++i) {
21: PetscMalloc1(stag->nRanks[i],&l[i]);
22: for (j=0; j<stag->nRanks[i]; ++j) l[i][j] = stag->l[i][j];
23: N[i] = stag->N[i];
24: }
26: /* dof */
27: dof = c < 0 ? -c : 1;
29: /* Determine/adjust sizes */
30: switch (loc) {
31: case DMSTAG_ELEMENT:
32: break;
33: case DMSTAG_LEFT:
34: case DMSTAG_RIGHT:
35: if (dim<1) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Incompatible dim (%d) and loc(%s) combination",dim,DMStagStencilLocations[loc]);
36: l[0][stag->nRanks[0]-1] += 1; /* extra vertex in direction 0 on last rank in dimension 0 */
37: N[0] += 1;
38: break;
39: case DMSTAG_UP:
40: case DMSTAG_DOWN:
41: if (dim < 2) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Incompatible dim (%d) and loc(%s) combination",dim,DMStagStencilLocations[loc]);
42: l[1][stag->nRanks[1]-1] += 1; /* extra vertex in direction 1 on last rank in dimension 1 */
43: N[1] += 1;
44: break;
45: case DMSTAG_BACK:
46: case DMSTAG_FRONT:
47: if (dim < 3) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Incompatible dim (%d) and loc(%s) combination",dim,DMStagStencilLocations[loc]);
48: l[2][stag->nRanks[2]-1] += 1; /* extra vertex in direction 2 on last rank in dimension 2 */
49: N[2] += 1;
50: break;
51: case DMSTAG_DOWN_LEFT :
52: case DMSTAG_DOWN_RIGHT :
53: case DMSTAG_UP_LEFT :
54: case DMSTAG_UP_RIGHT :
55: if (dim < 2) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Incompatible dim (%d) and loc(%s) combination",dim,DMStagStencilLocations[loc]);
56: for (i=0; i<2; ++i) { /* extra vertex in direction i on last rank in dimension i = 0,1 */
57: l[i][stag->nRanks[i]-1] += 1;
58: N[i] += 1;
59: }
60: break;
61: case DMSTAG_BACK_LEFT:
62: case DMSTAG_BACK_RIGHT:
63: case DMSTAG_FRONT_LEFT:
64: case DMSTAG_FRONT_RIGHT:
65: if (dim < 3) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Incompatible dim (%d) and loc(%s) combination",dim,DMStagStencilLocations[loc]);
66: for (i=0; i<3; i+=2) { /* extra vertex in direction i on last rank in dimension i = 0,2 */
67: l[i][stag->nRanks[i]-1] += 1;
68: N[i] += 1;
69: }
70: break;
71: case DMSTAG_BACK_DOWN:
72: case DMSTAG_BACK_UP:
73: case DMSTAG_FRONT_DOWN:
74: case DMSTAG_FRONT_UP:
75: if (dim < 3) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Incompatible dim (%d) and loc(%s) combination",dim,DMStagStencilLocations[loc]);
76: for (i=1; i<3; ++i) { /* extra vertex in direction i on last rank in dimension i = 1,2 */
77: l[i][stag->nRanks[i]-1] += 1;
78: N[i] += 1;
79: }
80: break;
81: case DMSTAG_BACK_DOWN_LEFT:
82: case DMSTAG_BACK_DOWN_RIGHT:
83: case DMSTAG_BACK_UP_LEFT:
84: case DMSTAG_BACK_UP_RIGHT:
85: case DMSTAG_FRONT_DOWN_LEFT:
86: case DMSTAG_FRONT_DOWN_RIGHT:
87: case DMSTAG_FRONT_UP_LEFT:
88: case DMSTAG_FRONT_UP_RIGHT:
89: if (dim < 3) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Incompatible dim (%d) and loc(%s) combination",dim,DMStagStencilLocations[loc]);
90: for (i=0; i<3; ++i) { /* extra vertex in direction i on last rank in dimension i = 0,1,2 */
91: l[i][stag->nRanks[i]-1] += 1;
92: N[i] += 1;
93: }
94: break;
95: break;
96: default : SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented for location %s",DMStagStencilLocations[loc]);
97: }
99: /* Use the same stencil type */
100: switch (stag->stencilType) {
101: case DMSTAG_STENCIL_STAR: stencilType = DMDA_STENCIL_STAR; stencilWidth = stag->stencilWidth; break;
102: case DMSTAG_STENCIL_BOX : stencilType = DMDA_STENCIL_BOX ; stencilWidth = stag->stencilWidth; break;
103: default: SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unsupported Stencil Type %d",stag->stencilType);
104: }
106: /* Create DMDA, using same boundary type */
107: switch (dim) {
108: case 1:
109: DMDACreate1d(PetscObjectComm((PetscObject)dm),stag->boundaryType[0],N[0],dof,stencilWidth,l[0],dmda);
110: break;
111: case 2:
112: DMDACreate2d(PetscObjectComm((PetscObject)dm),stag->boundaryType[0],stag->boundaryType[1],stencilType,N[0],N[1],stag->nRanks[0],stag->nRanks[1],dof,stencilWidth,l[0],l[1],dmda);
113: break;
114: case 3:
115: DMDACreate3d(PetscObjectComm((PetscObject)dm),stag->boundaryType[0],stag->boundaryType[1],stag->boundaryType[2],stencilType,N[0],N[1],N[2],stag->nRanks[0],stag->nRanks[1],stag->nRanks[2],dof,stencilWidth,l[0],l[1],l[2],dmda);
116: break;
117: default: SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not implemented for dim %d",dim);
118: }
119: for (i=0; i<dim; ++i) {
120: PetscFree(l[i]);
121: }
122: return(0);
123: }
125: /*
126: Helper function to get the number of extra points in a DMDA representation for a given canonical location.
127: */
128: static PetscErrorCode DMStagDMDAGetExtraPoints(DM dm,DMStagStencilLocation locCanonical,PetscInt *extraPoint)
129: {
131: PetscInt dim,d,nExtra[DMSTAG_MAX_DIM];
135: DMGetDimension(dm,&dim);
136: if (dim > DMSTAG_MAX_DIM) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented for %D dimensions",dim);
137: DMStagGetCorners(dm,NULL,NULL,NULL,NULL,NULL,NULL,&nExtra[0],&nExtra[1],&nExtra[2]);
138: for (d=0; d<dim; ++d) extraPoint[d] = 0;
139: switch (locCanonical) {
140: case DMSTAG_ELEMENT:
141: break; /* no extra points */
142: case DMSTAG_LEFT:
143: extraPoint[0] = nExtra[0]; break; /* only extra point in x */
144: case DMSTAG_DOWN:
145: extraPoint[1] = nExtra[1]; break; /* only extra point in y */
146: case DMSTAG_BACK:
147: extraPoint[2] = nExtra[2]; break; /* only extra point in z */
148: case DMSTAG_DOWN_LEFT:
149: extraPoint[0] = nExtra[0]; extraPoint[1] = nExtra[1]; break; /* extra point in both x and y */
150: case DMSTAG_BACK_LEFT:
151: extraPoint[0] = nExtra[0]; extraPoint[2] = nExtra[2]; break; /* extra point in both x and z */
152: case DMSTAG_BACK_DOWN:
153: extraPoint[1] = nExtra[1]; extraPoint[2] = nExtra[2]; break; /* extra point in both y and z */
154: case DMSTAG_BACK_DOWN_LEFT:
155: extraPoint[0] = nExtra[0]; extraPoint[1] = nExtra[1]; extraPoint[2] = nExtra[2]; break; /* extra points in x,y,z */
156: default : SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Not implemented for location (perhaps not canonical) %s",DMStagStencilLocations[locCanonical]);
157: }
158: return(0);
159: }
161: /*
162: Function much like DMStagMigrateVec(), but which accepts an additional position argument to disambiguate which
163: type of DMDA to migrate to.
164: */
166: static PetscErrorCode DMStagMigrateVecDMDA(DM dm,Vec vec,DMStagStencilLocation loc,PetscInt c,DM dmTo,Vec vecTo)
167: {
169: PetscInt i,j,k,d,dim,dof,dofToMax,start[DMSTAG_MAX_DIM],n[DMSTAG_MAX_DIM],extraPoint[DMSTAG_MAX_DIM];
170: Vec vecLocal;
177: DMGetDimension(dm,&dim);
178: DMDAGetDof(dmTo,&dofToMax);
179: if (-c > dofToMax) SETERRQ1(PetscObjectComm((PetscObject)dmTo),PETSC_ERR_ARG_OUTOFRANGE,"Invalid negative component value. Must be >= -%D",dofToMax);
180: DMStagGetCorners(dm,&start[0],&start[1],&start[2],&n[0],&n[1],&n[2],NULL,NULL,NULL);
181: DMStagDMDAGetExtraPoints(dm,loc,extraPoint);
182: DMStagGetLocationDOF(dm,loc,&dof);
183: DMGetLocalVector(dm,&vecLocal);
184: DMGlobalToLocalBegin(dm,vec,INSERT_VALUES,vecLocal);
185: DMGlobalToLocalEnd(dm,vec,INSERT_VALUES,vecLocal);
186: if (dim == 1) {
187: PetscScalar **arrTo;
188: DMDAVecGetArrayDOF(dmTo,vecTo,&arrTo);
189: if (c < 0) {
190: const PetscInt dofTo = -c;
191: for (i=start[0]; i<start[0] + n[0] + extraPoint[0]; ++i) {
192: for (d=0; d<PetscMin(dof,dofTo); ++d) {
193: DMStagStencil pos;
194: pos.i = i; pos.loc = loc; pos.c = d;
195: DMStagVecGetValuesStencil(dm,vecLocal,1,&pos,&arrTo[i][d]);
196: }
197: for (;d<dofTo; ++d) {
198: arrTo[i][d] = 0.0; /* Pad extra dof with zeros */
199: }
200: }
201: } else {
202: for (i=start[0]; i<start[0] + n[0] + extraPoint[0]; ++i) {
203: DMStagStencil pos;
204: pos.i = i; pos.loc = loc; pos.c = c;
205: DMStagVecGetValuesStencil(dm,vecLocal,1,&pos,&arrTo[i][0]);
206: }
207: }
208: DMDAVecRestoreArrayDOF(dmTo,vecTo,&arrTo);
209: } else if (dim == 2) {
210: PetscScalar ***arrTo;
211: DMDAVecGetArrayDOF(dmTo,vecTo,&arrTo);
212: if (c < 0) {
213: const PetscInt dofTo = -c;
214: for (j=start[1]; j<start[1] + n[1] + extraPoint[1]; ++j) {
215: for (i=start[0]; i<start[0] + n[0] + extraPoint[0]; ++i) {
216: for (d=0; d<PetscMin(dof,dofTo); ++d) {
217: DMStagStencil pos;
218: pos.i = i; pos.j = j; pos.loc = loc; pos.c = d;
219: DMStagVecGetValuesStencil(dm,vecLocal,1,&pos,&arrTo[j][i][d]);
220: }
221: for (;d<dofTo; ++d) {
222: arrTo[j][i][d] = 0.0; /* Pad extra dof with zeros */
223: }
224: }
225: }
226: } else {
227: for (j=start[1]; j<start[1] + n[1] + extraPoint[1]; ++j) {
228: for (i=start[0]; i<start[0] + n[0] + extraPoint[0]; ++i) {
229: DMStagStencil pos;
230: pos.i = i; pos.j = j; pos.loc = loc; pos.c = c;
231: DMStagVecGetValuesStencil(dm,vecLocal,1,&pos,&arrTo[j][i][0]);
232: }
233: }
234: }
235: DMDAVecRestoreArrayDOF(dmTo,vecTo,&arrTo);
236: } else if (dim == 3) {
237: PetscScalar ****arrTo;
238: DMDAVecGetArrayDOF(dmTo,vecTo,&arrTo);
239: if (c < 0) {
240: const PetscInt dofTo = -c;
241: for (k=start[2]; k<start[2] + n[2] + extraPoint[2]; ++k) {
242: for (j=start[1]; j<start[1] + n[1] + extraPoint[1]; ++j) {
243: for (i=start[0]; i<start[0] + n[0] + extraPoint[0]; ++i) {
244: for (d=0; d<PetscMin(dof,dofTo); ++d) {
245: DMStagStencil pos;
246: pos.i = i; pos.j = j; pos.k = k; pos.loc = loc; pos.c = d;
247: DMStagVecGetValuesStencil(dm,vecLocal,1,&pos,&arrTo[k][j][i][d]);
248: }
249: for (;d<dofTo; ++d) {
250: arrTo[k][j][i][d] = 0.0; /* Pad extra dof with zeros */
251: }
252: }
253: }
254: }
255: } else {
256: for (k=start[2]; k<start[2] + n[2] + extraPoint[2]; ++k) {
257: for (j=start[1]; j<start[1] + n[1] + extraPoint[1]; ++j) {
258: for (i=start[0]; i<start[0] + n[0] + extraPoint[0]; ++i) {
259: DMStagStencil pos;
260: pos.i = i; pos.j = j; pos.k = k; pos.loc = loc; pos.c = c;
261: DMStagVecGetValuesStencil(dm,vecLocal,1,&pos,&arrTo[k][j][i][0]);
262: }
263: }
264: }
265: }
266: DMDAVecRestoreArrayDOF(dmTo,vecTo,&arrTo);
267: } else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unsupported dimension %d",dim);
268: DMRestoreLocalVector(dm,&vecLocal);
269: return(0);
270: }
272: /* Transfer coordinates from a DMStag to a DMDA, specifying which location */
273: static PetscErrorCode DMStagTransferCoordinatesToDMDA(DM dmstag,DMStagStencilLocation loc,DM dmda)
274: {
276: PetscInt dim,start[DMSTAG_MAX_DIM],n[DMSTAG_MAX_DIM],extraPoint[DMSTAG_MAX_DIM],d;
277: DM dmstagCoord,dmdaCoord;
278: DMType dmstagCoordType;
279: Vec stagCoord,daCoord;
280: PetscBool daCoordIsStag,daCoordIsProduct;
285: DMGetDimension(dmstag,&dim);
286: DMGetCoordinateDM(dmstag,&dmstagCoord);
287: DMGetCoordinatesLocal(dmstag,&stagCoord); /* Note local */
288: DMGetCoordinateDM(dmda,&dmdaCoord);
289: daCoord = NULL;
290: DMGetCoordinates(dmda,&daCoord);
291: if (!daCoord) {
292: DMCreateGlobalVector(dmdaCoord,&daCoord);
293: DMSetCoordinates(dmda,daCoord);
294: VecDestroy(&daCoord);
295: DMGetCoordinates(dmda,&daCoord);
296: }
297: DMGetType(dmstagCoord,&dmstagCoordType);
298: PetscStrcmp(dmstagCoordType,DMSTAG,&daCoordIsStag);
299: PetscStrcmp(dmstagCoordType,DMPRODUCT,&daCoordIsProduct);
300: DMStagGetCorners(dmstag,&start[0],&start[1],&start[2],&n[0],&n[1],&n[2],NULL,NULL,NULL);
301: DMStagDMDAGetExtraPoints(dmstag,loc,extraPoint);
302: if (dim == 1) {
303: PetscInt ex;
304: PetscScalar **cArrDa;
305: DMDAVecGetArrayDOF(dmdaCoord,daCoord,&cArrDa);
306: if (daCoordIsStag) {
307: PetscInt slot;
308: PetscScalar **cArrStag;
309: DMStagGetLocationSlot(dmstagCoord,loc,0,&slot);
310: DMStagVecGetArrayRead(dmstagCoord,stagCoord,&cArrStag);
311: for (ex=start[0]; ex<start[0] + n[0] + extraPoint[0]; ++ex) {
312: cArrDa[ex][0] = cArrStag[ex][slot];
313: }
314: DMStagVecRestoreArrayRead(dmstagCoord,stagCoord,&cArrStag);
315: } else if (daCoordIsProduct) {
316: PetscScalar **cArrX;
317: DMStagGetProductCoordinateArraysRead(dmstag,&cArrX,NULL,NULL);
318: for (ex=start[0]; ex<start[0] + n[0] + extraPoint[0]; ++ex) {
319: cArrDa[ex][0] = cArrX[ex][0];
320: }
321: DMStagRestoreProductCoordinateArraysRead(dmstag,&cArrX,NULL,NULL);
322: } else SETERRQ(PetscObjectComm((PetscObject)dmstag),PETSC_ERR_SUP,"Stag to DA coordinate transfer only supported for DMStag coordinate DM of type DMstag or DMProduct");
323: DMDAVecRestoreArrayDOF(dmdaCoord,daCoord,&cArrDa);
324: } else if (dim == 2) {
325: PetscInt ex,ey;
326: PetscScalar ***cArrDa;
327: DMDAVecGetArrayDOF(dmdaCoord,daCoord,&cArrDa);
328: if (daCoordIsStag) {
329: PetscInt slot;
330: PetscScalar ***cArrStag;
331: DMStagGetLocationSlot(dmstagCoord,loc,0,&slot);
332: DMStagVecGetArrayRead(dmstagCoord,stagCoord,&cArrStag);
333: for (ey=start[1]; ey<start[1] + n[1] + extraPoint[1]; ++ey) {
334: for (ex=start[0]; ex<start[0] + n[0] + extraPoint[0]; ++ex) {
335: for (d=0; d<2; ++d) {
336: cArrDa[ey][ex][d] = cArrStag[ey][ex][slot+d];
337: }
338: }
339: }
340: DMStagVecRestoreArrayRead(dmstagCoord,stagCoord,&cArrStag);
341: } else if (daCoordIsProduct) {
342: PetscScalar **cArrX,**cArrY;
343: DMStagGetProductCoordinateArraysRead(dmstag,&cArrX,&cArrY,NULL);
344: for (ey=start[1]; ey<start[1] + n[1] + extraPoint[1]; ++ey) {
345: for (ex=start[0]; ex<start[0] + n[0] + extraPoint[0]; ++ex) {
346: cArrDa[ey][ex][0] = cArrX[ex][0];
347: cArrDa[ey][ex][1] = cArrY[ey][0];
348: }
349: }
350: DMStagRestoreProductCoordinateArraysRead(dmstag,&cArrX,&cArrY,NULL);
351: } else SETERRQ(PetscObjectComm((PetscObject)dmstag),PETSC_ERR_SUP,"Stag to DA coordinate transfer only supported for DMStag coordinate DM of type DMstag or DMProduct");
352: DMDAVecRestoreArrayDOF(dmdaCoord,daCoord,&cArrDa);
353: } else if (dim == 3) {
354: PetscInt ex,ey,ez;
355: PetscScalar ****cArrDa;
356: DMDAVecGetArrayDOF(dmdaCoord,daCoord,&cArrDa);
357: if (daCoordIsStag) {
358: PetscInt slot;
359: PetscScalar ****cArrStag;
360: DMStagGetLocationSlot(dmstagCoord,loc,0,&slot);
361: DMStagVecGetArrayRead(dmstagCoord,stagCoord,&cArrStag);
362: for (ez=start[2]; ez<start[2] + n[2] + extraPoint[2]; ++ez) {
363: for (ey=start[1]; ey<start[1] + n[1] + extraPoint[1]; ++ey) {
364: for (ex=start[0]; ex<start[0] + n[0] + extraPoint[0]; ++ex) {
365: for (d=0; d<3; ++d) {
366: cArrDa[ez][ey][ex][d] = cArrStag[ez][ey][ex][slot+d];
367: }
368: }
369: }
370: }
371: DMStagVecRestoreArrayRead(dmstagCoord,stagCoord,&cArrStag);
372: } else if (daCoordIsProduct) {
373: PetscScalar **cArrX,**cArrY,**cArrZ;
374: DMStagGetProductCoordinateArraysRead(dmstag,&cArrX,&cArrY,&cArrZ);
375: for (ez=start[2]; ez<start[2] + n[2] + extraPoint[2]; ++ez) {
376: for (ey=start[1]; ey<start[1] + n[1] + extraPoint[1]; ++ey) {
377: for (ex=start[0]; ex<start[0] + n[0] + extraPoint[0]; ++ex) {
378: cArrDa[ez][ey][ex][0] = cArrX[ex][0];
379: cArrDa[ez][ey][ex][1] = cArrY[ey][0];
380: cArrDa[ez][ey][ex][2] = cArrZ[ez][0];
381: }
382: }
383: }
384: DMStagRestoreProductCoordinateArraysRead(dmstag,&cArrX,&cArrY,&cArrZ);
385: } else SETERRQ(PetscObjectComm((PetscObject)dmstag),PETSC_ERR_SUP,"Stag to DA coordinate transfer only supported for DMStag coordinate DM of type DMstag or DMProduct");
386: DMDAVecRestoreArrayDOF(dmdaCoord,daCoord,&cArrDa);
387: } else SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Unsupported dimension %d",dim);
388: return(0);
389: }
391: /*
392: Convert to a location value with only BACK, DOWN, LEFT, and ELEMENT involved (makes looping easier)
393: */
394: static PetscErrorCode DMStagStencilLocationCanonicalize(DMStagStencilLocation loc,DMStagStencilLocation *locCanonical)
395: {
397: switch (loc) {
398: case DMSTAG_ELEMENT:
399: *locCanonical = DMSTAG_ELEMENT;
400: break;
401: case DMSTAG_LEFT:
402: case DMSTAG_RIGHT:
403: *locCanonical = DMSTAG_LEFT;
404: break;
405: case DMSTAG_DOWN:
406: case DMSTAG_UP:
407: *locCanonical = DMSTAG_DOWN;
408: break;
409: case DMSTAG_BACK:
410: case DMSTAG_FRONT:
411: *locCanonical = DMSTAG_BACK;
412: break;
413: case DMSTAG_DOWN_LEFT :
414: case DMSTAG_DOWN_RIGHT :
415: case DMSTAG_UP_LEFT :
416: case DMSTAG_UP_RIGHT :
417: *locCanonical = DMSTAG_DOWN_LEFT;
418: break;
419: case DMSTAG_BACK_LEFT:
420: case DMSTAG_BACK_RIGHT:
421: case DMSTAG_FRONT_LEFT:
422: case DMSTAG_FRONT_RIGHT:
423: *locCanonical = DMSTAG_BACK_LEFT;
424: break;
425: case DMSTAG_BACK_DOWN:
426: case DMSTAG_BACK_UP:
427: case DMSTAG_FRONT_DOWN:
428: case DMSTAG_FRONT_UP:
429: *locCanonical = DMSTAG_BACK_DOWN;
430: break;
431: case DMSTAG_BACK_DOWN_LEFT:
432: case DMSTAG_BACK_DOWN_RIGHT:
433: case DMSTAG_BACK_UP_LEFT:
434: case DMSTAG_BACK_UP_RIGHT:
435: case DMSTAG_FRONT_DOWN_LEFT:
436: case DMSTAG_FRONT_DOWN_RIGHT:
437: case DMSTAG_FRONT_UP_LEFT:
438: case DMSTAG_FRONT_UP_RIGHT:
439: *locCanonical = DMSTAG_BACK_DOWN_LEFT;
440: break;
441: default :
442: *locCanonical = DMSTAG_NULL_LOCATION;
443: break;
444: }
445: return(0);
446: }
448: /*@C
449: DMStagVecSplitToDMDA - create a DMDA and Vec from a DMStag and Vec
451: Logically Collective
453: High-level helper function which accepts a DMStag, a global vector, and location/dof,
454: and generates a corresponding DMDA and Vec.
456: Input Parameters:
457: + dm - the DMStag object
458: . vec- Vec object associated with dm
459: . loc - which subgrid to extract (see DMStagStencilLocation)
460: - c - which component to extract (see note below)
462: Output Parameters:
463: + pda - the new DMDA
464: - pdavec - the new Vec
466: Notes:
467: If a c value of -k is provided, the first k dof for that position are extracted,
468: padding with zero values if needbe. If a non-negative value is provided, a single
469: dof is extracted.
471: The caller is responsible for destroying the created DMDA and Vec.
473: Level: advanced
475: .seealso: DMSTAG, DMDA, DMStagMigrateVec(), DMStagCreateCompatibleDMStag()
476: @*/
477: PetscErrorCode DMStagVecSplitToDMDA(DM dm,Vec vec,DMStagStencilLocation loc,PetscInt c,DM *pda,Vec *pdavec)
478: {
479: PetscErrorCode ierr;
480: PetscInt dim,locdof;
481: DM da,coordDM;
482: Vec davec;
483: DMStagStencilLocation locCanonical;
488: DMGetDimension(dm,&dim);
489: DMStagGetLocationDOF(dm,loc,&locdof);
490: if (c >= locdof) SETERRQ3(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Location %s has %D dof, but component %D requested\n",DMStagStencilLocations[loc],locdof,c);
491: DMStagStencilLocationCanonicalize(loc,&locCanonical);
492: DMStagCreateCompatibleDMDA(dm,locCanonical,c,pda);
493: da = *pda;
494: DMSetUp(*pda);
495: DMGetCoordinateDM(dm,&coordDM);
496: if (coordDM) {
497: DMStagTransferCoordinatesToDMDA(dm,locCanonical,da);
498: }
499: DMCreateGlobalVector(da,pdavec);
500: davec = *pdavec;
501: DMStagMigrateVecDMDA(dm,vec,locCanonical,c,da,davec);
502: return(0);
503: }