Actual source code: plexinterpolate.c
petsc-3.14.6 2021-03-30
1: #include <petsc/private/dmpleximpl.h>
2: #include <petsc/private/hashmapi.h>
3: #include <petsc/private/hashmapij.h>
5: const char * const DMPlexInterpolatedFlags[] = {"none", "partial", "mixed", "full", "DMPlexInterpolatedFlag", "DMPLEX_INTERPOLATED_", NULL};
7: /* HashIJKL */
9: #include <petsc/private/hashmap.h>
11: typedef struct _PetscHashIJKLKey { PetscInt i, j, k, l; } PetscHashIJKLKey;
13: #define PetscHashIJKLKeyHash(key) \
14: PetscHashCombine(PetscHashCombine(PetscHashInt((key).i),PetscHashInt((key).j)), \
15: PetscHashCombine(PetscHashInt((key).k),PetscHashInt((key).l)))
17: #define PetscHashIJKLKeyEqual(k1,k2) \
18: (((k1).i==(k2).i) ? ((k1).j==(k2).j) ? ((k1).k==(k2).k) ? ((k1).l==(k2).l) : 0 : 0 : 0)
20: PETSC_HASH_MAP(HashIJKL, PetscHashIJKLKey, PetscInt, PetscHashIJKLKeyHash, PetscHashIJKLKeyEqual, -1)
22: static PetscSFNode _PetscInvalidSFNode = {-1, -1};
24: typedef struct _PetscHashIJKLRemoteKey { PetscSFNode i, j, k, l; } PetscHashIJKLRemoteKey;
26: #define PetscHashIJKLRemoteKeyHash(key) \
27: PetscHashCombine(PetscHashCombine(PetscHashInt((key).i.rank + (key).i.index),PetscHashInt((key).j.rank + (key).j.index)), \
28: PetscHashCombine(PetscHashInt((key).k.rank + (key).k.index),PetscHashInt((key).l.rank + (key).l.index)))
30: #define PetscHashIJKLRemoteKeyEqual(k1,k2) \
31: (((k1).i.rank==(k2).i.rank) ? ((k1).i.index==(k2).i.index) ? ((k1).j.rank==(k2).j.rank) ? ((k1).j.index==(k2).j.index) ? ((k1).k.rank==(k2).k.rank) ? ((k1).k.index==(k2).k.index) ? ((k1).l.rank==(k2).l.rank) ? ((k1).l.index==(k2).l.index) : 0 : 0 : 0 : 0 : 0 : 0 : 0)
33: PETSC_HASH_MAP(HashIJKLRemote, PetscHashIJKLRemoteKey, PetscSFNode, PetscHashIJKLRemoteKeyHash, PetscHashIJKLRemoteKeyEqual, _PetscInvalidSFNode)
35: static PetscErrorCode PetscSortSFNode(PetscInt n, PetscSFNode A[])
36: {
37: PetscInt i;
40: for (i = 1; i < n; ++i) {
41: PetscSFNode x = A[i];
42: PetscInt j;
44: for (j = i-1; j >= 0; --j) {
45: if ((A[j].rank > x.rank) || (A[j].rank == x.rank && A[j].index > x.index)) break;
46: A[j+1] = A[j];
47: }
48: A[j+1] = x;
49: }
50: return(0);
51: }
53: /*
54: DMPlexGetRawFaces_Internal - Gets groups of vertices that correspond to faces for the given cone
55: */
56: PetscErrorCode DMPlexGetRawFaces_Internal(DM dm, DMPolytopeType ct, const PetscInt cone[], PetscInt *numFaces, const DMPolytopeType *faceTypes[], const PetscInt *faceSizes[], const PetscInt *faces[])
57: {
58: DMPolytopeType *typesTmp;
59: PetscInt *sizesTmp, *facesTmp;
60: PetscInt maxConeSize, maxSupportSize;
61: PetscErrorCode ierr;
66: DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);
67: if (faceTypes) {DMGetWorkArray(dm, PetscMax(maxConeSize, maxSupportSize), MPIU_INT, &typesTmp);}
68: if (faceSizes) {DMGetWorkArray(dm, PetscMax(maxConeSize, maxSupportSize), MPIU_INT, &sizesTmp);}
69: if (faces) {DMGetWorkArray(dm, PetscSqr(PetscMax(maxConeSize, maxSupportSize)), MPIU_INT, &facesTmp);}
70: switch (ct) {
71: case DM_POLYTOPE_POINT:
72: if (numFaces) *numFaces = 0;
73: break;
74: case DM_POLYTOPE_SEGMENT:
75: if (numFaces) *numFaces = 2;
76: if (faceTypes) {
77: typesTmp[0] = DM_POLYTOPE_POINT; typesTmp[1] = DM_POLYTOPE_POINT;
78: *faceTypes = typesTmp;
79: }
80: if (faceSizes) {
81: sizesTmp[0] = 1; sizesTmp[1] = 1;
82: *faceSizes = sizesTmp;
83: }
84: if (faces) {
85: facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
86: *faces = facesTmp;
87: }
88: break;
89: case DM_POLYTOPE_POINT_PRISM_TENSOR:
90: if (numFaces) *numFaces = 2;
91: if (faceTypes) {
92: typesTmp[0] = DM_POLYTOPE_POINT; typesTmp[1] = DM_POLYTOPE_POINT;
93: *faceTypes = typesTmp;
94: }
95: if (faceSizes) {
96: sizesTmp[0] = 1; sizesTmp[1] = 1;
97: *faceSizes = sizesTmp;
98: }
99: if (faces) {
100: facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
101: *faces = facesTmp;
102: }
103: break;
104: case DM_POLYTOPE_TRIANGLE:
105: if (numFaces) *numFaces = 3;
106: if (faceTypes) {
107: typesTmp[0] = DM_POLYTOPE_SEGMENT; typesTmp[1] = DM_POLYTOPE_SEGMENT; typesTmp[2] = DM_POLYTOPE_SEGMENT;
108: *faceTypes = typesTmp;
109: }
110: if (faceSizes) {
111: sizesTmp[0] = 2; sizesTmp[1] = 2; sizesTmp[2] = 2;
112: *faceSizes = sizesTmp;
113: }
114: if (faces) {
115: facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
116: facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
117: facesTmp[4] = cone[2]; facesTmp[5] = cone[0];
118: *faces = facesTmp;
119: }
120: break;
121: case DM_POLYTOPE_QUADRILATERAL:
122: /* Vertices follow right hand rule */
123: if (numFaces) *numFaces = 4;
124: if (faceTypes) {
125: typesTmp[0] = DM_POLYTOPE_SEGMENT; typesTmp[1] = DM_POLYTOPE_SEGMENT; typesTmp[2] = DM_POLYTOPE_SEGMENT; typesTmp[3] = DM_POLYTOPE_SEGMENT;
126: *faceTypes = typesTmp;
127: }
128: if (faceSizes) {
129: sizesTmp[0] = 2; sizesTmp[1] = 2; sizesTmp[2] = 2; sizesTmp[3] = 2;
130: *faceSizes = sizesTmp;
131: }
132: if (faces) {
133: facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
134: facesTmp[2] = cone[1]; facesTmp[3] = cone[2];
135: facesTmp[4] = cone[2]; facesTmp[5] = cone[3];
136: facesTmp[6] = cone[3]; facesTmp[7] = cone[0];
137: *faces = facesTmp;
138: }
139: break;
140: case DM_POLYTOPE_SEG_PRISM_TENSOR:
141: if (numFaces) *numFaces = 4;
142: if (faceTypes) {
143: typesTmp[0] = DM_POLYTOPE_SEGMENT; typesTmp[1] = DM_POLYTOPE_SEGMENT; typesTmp[2] = DM_POLYTOPE_POINT_PRISM_TENSOR; typesTmp[3] = DM_POLYTOPE_POINT_PRISM_TENSOR;
144: *faceTypes = typesTmp;
145: }
146: if (faceSizes) {
147: sizesTmp[0] = 2; sizesTmp[1] = 2; sizesTmp[2] = 2; sizesTmp[3] = 2;
148: *faceSizes = sizesTmp;
149: }
150: if (faces) {
151: facesTmp[0] = cone[0]; facesTmp[1] = cone[1];
152: facesTmp[2] = cone[2]; facesTmp[3] = cone[3];
153: facesTmp[4] = cone[0]; facesTmp[5] = cone[2];
154: facesTmp[6] = cone[1]; facesTmp[7] = cone[3];
155: *faces = facesTmp;
156: }
157: break;
158: case DM_POLYTOPE_TETRAHEDRON:
159: /* Vertices of first face follow right hand rule and normal points away from last vertex */
160: if (numFaces) *numFaces = 4;
161: if (faceTypes) {
162: typesTmp[0] = DM_POLYTOPE_TRIANGLE; typesTmp[1] = DM_POLYTOPE_TRIANGLE; typesTmp[2] = DM_POLYTOPE_TRIANGLE; typesTmp[3] = DM_POLYTOPE_TRIANGLE;
163: *faceTypes = typesTmp;
164: }
165: if (faceSizes) {
166: sizesTmp[0] = 3; sizesTmp[1] = 3; sizesTmp[2] = 3; sizesTmp[3] = 3;
167: *faceSizes = sizesTmp;
168: }
169: if (faces) {
170: facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2];
171: facesTmp[3] = cone[0]; facesTmp[4] = cone[3]; facesTmp[5] = cone[1];
172: facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[3];
173: facesTmp[9] = cone[2]; facesTmp[10] = cone[1]; facesTmp[11] = cone[3];
174: *faces = facesTmp;
175: }
176: break;
177: case DM_POLYTOPE_HEXAHEDRON:
178: /* 7--------6
179: /| /|
180: / | / |
181: 4--------5 |
182: | | | |
183: | | | |
184: | 1--------2
185: | / | /
186: |/ |/
187: 0--------3
188: */
189: if (numFaces) *numFaces = 6;
190: if (faceTypes) {
191: typesTmp[0] = DM_POLYTOPE_QUADRILATERAL; typesTmp[1] = DM_POLYTOPE_QUADRILATERAL; typesTmp[2] = DM_POLYTOPE_QUADRILATERAL;
192: typesTmp[3] = DM_POLYTOPE_QUADRILATERAL; typesTmp[4] = DM_POLYTOPE_QUADRILATERAL; typesTmp[5] = DM_POLYTOPE_QUADRILATERAL;
193: *faceTypes = typesTmp;
194: }
195: if (faceSizes) {
196: sizesTmp[0] = 4; sizesTmp[1] = 4; sizesTmp[2] = 4; sizesTmp[3] = 4; sizesTmp[4] = 4; sizesTmp[5] = 4;
197: *faceSizes = sizesTmp;
198: }
199: if (faces) {
200: facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; /* Bottom */
201: facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; /* Top */
202: facesTmp[8] = cone[0]; facesTmp[9] = cone[3]; facesTmp[10] = cone[5]; facesTmp[11] = cone[4]; /* Front */
203: facesTmp[12] = cone[2]; facesTmp[13] = cone[1]; facesTmp[14] = cone[7]; facesTmp[15] = cone[6]; /* Back */
204: facesTmp[16] = cone[3]; facesTmp[17] = cone[2]; facesTmp[18] = cone[6]; facesTmp[19] = cone[5]; /* Right */
205: facesTmp[20] = cone[0]; facesTmp[21] = cone[4]; facesTmp[22] = cone[7]; facesTmp[23] = cone[1]; /* Left */
206: *faces = facesTmp;
207: }
208: break;
209: case DM_POLYTOPE_TRI_PRISM:
210: if (numFaces) *numFaces = 5;
211: if (faceTypes) {
212: typesTmp[0] = DM_POLYTOPE_TRIANGLE; typesTmp[1] = DM_POLYTOPE_TRIANGLE;
213: typesTmp[2] = DM_POLYTOPE_QUADRILATERAL; typesTmp[3] = DM_POLYTOPE_QUADRILATERAL; typesTmp[4] = DM_POLYTOPE_QUADRILATERAL;
214: *faceTypes = typesTmp;
215: }
216: if (faceSizes) {
217: sizesTmp[0] = 3; sizesTmp[1] = 3;
218: sizesTmp[2] = 4; sizesTmp[3] = 4; sizesTmp[4] = 4;
219: *faceSizes = sizesTmp;
220: }
221: if (faces) {
222: facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; /* Bottom */
223: facesTmp[3] = cone[3]; facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; /* Top */
224: facesTmp[6] = cone[0]; facesTmp[7] = cone[2]; facesTmp[8] = cone[4]; facesTmp[9] = cone[3]; /* Back left */
225: facesTmp[10] = cone[2]; facesTmp[11] = cone[1]; facesTmp[12] = cone[5]; facesTmp[13] = cone[4]; /* Front */
226: facesTmp[14] = cone[1]; facesTmp[15] = cone[0]; facesTmp[16] = cone[3]; facesTmp[17] = cone[5]; /* Back right */
227: *faces = facesTmp;
228: }
229: break;
230: case DM_POLYTOPE_TRI_PRISM_TENSOR:
231: if (numFaces) *numFaces = 5;
232: if (faceTypes) {
233: typesTmp[0] = DM_POLYTOPE_TRIANGLE; typesTmp[1] = DM_POLYTOPE_TRIANGLE;
234: typesTmp[2] = DM_POLYTOPE_SEG_PRISM_TENSOR; typesTmp[3] = DM_POLYTOPE_SEG_PRISM_TENSOR; typesTmp[4] = DM_POLYTOPE_SEG_PRISM_TENSOR;
235: *faceTypes = typesTmp;
236: }
237: if (faceSizes) {
238: sizesTmp[0] = 3; sizesTmp[1] = 3;
239: sizesTmp[2] = 4; sizesTmp[3] = 4; sizesTmp[4] = 4;
240: *faceSizes = sizesTmp;
241: }
242: if (faces) {
243: facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; /* Bottom */
244: facesTmp[3] = cone[3]; facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; /* Top */
245: facesTmp[6] = cone[0]; facesTmp[7] = cone[1]; facesTmp[8] = cone[3]; facesTmp[9] = cone[4]; /* Back left */
246: facesTmp[10] = cone[1]; facesTmp[11] = cone[2]; facesTmp[12] = cone[4]; facesTmp[13] = cone[5]; /* Back right */
247: facesTmp[14] = cone[2]; facesTmp[15] = cone[0]; facesTmp[16] = cone[5]; facesTmp[17] = cone[3]; /* Front */
248: *faces = facesTmp;
249: }
250: break;
251: case DM_POLYTOPE_QUAD_PRISM_TENSOR:
252: /* 7--------6
253: /| /|
254: / | / |
255: 4--------5 |
256: | | | |
257: | | | |
258: | 3--------2
259: | / | /
260: |/ |/
261: 0--------1
262: */
263: if (numFaces) *numFaces = 6;
264: if (faceTypes) {
265: typesTmp[0] = DM_POLYTOPE_QUADRILATERAL; typesTmp[1] = DM_POLYTOPE_QUADRILATERAL; typesTmp[2] = DM_POLYTOPE_SEG_PRISM_TENSOR;
266: typesTmp[3] = DM_POLYTOPE_SEG_PRISM_TENSOR; typesTmp[4] = DM_POLYTOPE_SEG_PRISM_TENSOR; typesTmp[5] = DM_POLYTOPE_SEG_PRISM_TENSOR;
267: *faceTypes = typesTmp;
268: }
269: if (faceSizes) {
270: sizesTmp[0] = 4; sizesTmp[1] = 4; sizesTmp[2] = 4; sizesTmp[3] = 4; sizesTmp[4] = 4; sizesTmp[5] = 4;
271: *faceSizes = sizesTmp;
272: }
273: if (faces) {
274: facesTmp[0] = cone[0]; facesTmp[1] = cone[1]; facesTmp[2] = cone[2]; facesTmp[3] = cone[3]; /* Bottom */
275: facesTmp[4] = cone[4]; facesTmp[5] = cone[5]; facesTmp[6] = cone[6]; facesTmp[7] = cone[7]; /* Top */
276: facesTmp[8] = cone[0]; facesTmp[9] = cone[1]; facesTmp[10] = cone[4]; facesTmp[11] = cone[5]; /* Front */
277: facesTmp[12] = cone[1]; facesTmp[13] = cone[2]; facesTmp[14] = cone[5]; facesTmp[15] = cone[6]; /* Right */
278: facesTmp[16] = cone[2]; facesTmp[17] = cone[3]; facesTmp[18] = cone[6]; facesTmp[19] = cone[7]; /* Back */
279: facesTmp[20] = cone[3]; facesTmp[21] = cone[0]; facesTmp[22] = cone[7]; facesTmp[23] = cone[4]; /* Left */
280: *faces = facesTmp;
281: }
282: break;
283: default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No face description for cell type %s", DMPolytopeTypes[ct]);
284: }
285: return(0);
286: }
288: PetscErrorCode DMPlexRestoreRawFaces_Internal(DM dm, DMPolytopeType ct, const PetscInt cone[], PetscInt *numFaces, const DMPolytopeType *faceTypes[], const PetscInt *faceSizes[], const PetscInt *faces[])
289: {
290: PetscErrorCode ierr;
293: if (faceTypes) {DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faceTypes);}
294: if (faceSizes) {DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faceSizes);}
295: if (faces) {DMRestoreWorkArray(dm, 0, MPIU_INT, (void *) faces);}
296: return(0);
297: }
299: /* This interpolates faces for cells at some stratum */
300: static PetscErrorCode DMPlexInterpolateFaces_Internal(DM dm, PetscInt cellDepth, DM idm)
301: {
302: DMLabel ctLabel;
303: PetscHashIJKL faceTable;
304: PetscInt faceTypeNum[DM_NUM_POLYTOPES];
305: PetscInt depth, d, Np, cStart, cEnd, c, fStart, fEnd;
309: DMPlexGetDepth(dm, &depth);
310: PetscHashIJKLCreate(&faceTable);
311: PetscArrayzero(faceTypeNum, DM_NUM_POLYTOPES);
312: DMPlexGetDepthStratum(dm, cellDepth, &cStart, &cEnd);
313: /* Number new faces and save face vertices in hash table */
314: DMPlexGetDepthStratum(dm, depth > cellDepth ? cellDepth : 0, NULL, &fStart);
315: fEnd = fStart;
316: for (c = cStart; c < cEnd; ++c) {
317: const PetscInt *cone, *faceSizes, *faces;
318: const DMPolytopeType *faceTypes;
319: DMPolytopeType ct;
320: PetscInt numFaces, cf, foff = 0;
322: DMPlexGetCellType(dm, c, &ct);
323: DMPlexGetCone(dm, c, &cone);
324: DMPlexGetRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);
325: for (cf = 0; cf < numFaces; foff += faceSizes[cf], ++cf) {
326: const PetscInt faceSize = faceSizes[cf];
327: const DMPolytopeType faceType = faceTypes[cf];
328: const PetscInt *face = &faces[foff];
329: PetscHashIJKLKey key;
330: PetscHashIter iter;
331: PetscBool missing;
333: if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not support faces of size %D > 4", faceSize);
334: key.i = face[0];
335: key.j = faceSize > 1 ? face[1] : PETSC_MAX_INT;
336: key.k = faceSize > 2 ? face[2] : PETSC_MAX_INT;
337: key.l = faceSize > 3 ? face[3] : PETSC_MAX_INT;
338: PetscSortInt(faceSize, (PetscInt *) &key);
339: PetscHashIJKLPut(faceTable, key, &iter, &missing);
340: if (missing) {
341: PetscHashIJKLIterSet(faceTable, iter, fEnd++);
342: ++faceTypeNum[faceType];
343: }
344: }
345: DMPlexRestoreRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);
346: }
347: /* We need to number faces contiguously among types */
348: {
349: PetscInt faceTypeStart[DM_NUM_POLYTOPES], ct, numFT = 0;
351: for (ct = 0; ct < DM_NUM_POLYTOPES; ++ct) {if (faceTypeNum[ct]) ++numFT; faceTypeStart[ct] = 0;}
352: if (numFT > 1) {
353: PetscHashIJKLClear(faceTable);
354: faceTypeStart[0] = fStart;
355: for (ct = 1; ct < DM_NUM_POLYTOPES; ++ct) faceTypeStart[ct] = faceTypeStart[ct-1] + faceTypeNum[ct-1];
356: for (c = cStart; c < cEnd; ++c) {
357: const PetscInt *cone, *faceSizes, *faces;
358: const DMPolytopeType *faceTypes;
359: DMPolytopeType ct;
360: PetscInt numFaces, cf, foff = 0;
362: DMPlexGetCellType(dm, c, &ct);
363: DMPlexGetCone(dm, c, &cone);
364: DMPlexGetRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);
365: for (cf = 0; cf < numFaces; foff += faceSizes[cf], ++cf) {
366: const PetscInt faceSize = faceSizes[cf];
367: const DMPolytopeType faceType = faceTypes[cf];
368: const PetscInt *face = &faces[foff];
369: PetscHashIJKLKey key;
370: PetscHashIter iter;
371: PetscBool missing;
373: if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not support faces of size %D > 4", faceSize);
374: key.i = face[0];
375: key.j = faceSize > 1 ? face[1] : PETSC_MAX_INT;
376: key.k = faceSize > 2 ? face[2] : PETSC_MAX_INT;
377: key.l = faceSize > 3 ? face[3] : PETSC_MAX_INT;
378: PetscSortInt(faceSize, (PetscInt *) &key);
379: PetscHashIJKLPut(faceTable, key, &iter, &missing);
380: if (missing) {PetscHashIJKLIterSet(faceTable, iter, faceTypeStart[faceType]++);}
381: }
382: DMPlexRestoreRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);
383: }
384: for (ct = 1; ct < DM_NUM_POLYTOPES; ++ct) {
385: if (faceTypeStart[ct] != faceTypeStart[ct-1] + faceTypeNum[ct]) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent numbering for cell type %s, %D != %D + %D", DMPolytopeTypes[ct], faceTypeStart[ct], faceTypeStart[ct-1], faceTypeNum[ct]);
386: }
387: }
388: }
389: /* Add new points, always at the end of the numbering */
390: DMPlexGetChart(dm, NULL, &Np);
391: DMPlexSetChart(idm, 0, Np + (fEnd - fStart));
392: /* Set cone sizes */
393: /* Must create the celltype label here so that we do not automatically try to compute the types */
394: DMCreateLabel(idm, "celltype");
395: DMPlexGetCellTypeLabel(idm, &ctLabel);
396: for (d = 0; d <= depth; ++d) {
397: DMPolytopeType ct;
398: PetscInt coneSize, pStart, pEnd, p;
400: if (d == cellDepth) continue;
401: DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
402: for (p = pStart; p < pEnd; ++p) {
403: DMPlexGetConeSize(dm, p, &coneSize);
404: DMPlexSetConeSize(idm, p, coneSize);
405: DMPlexGetCellType(dm, p, &ct);
406: DMPlexSetCellType(idm, p, ct);
407: }
408: }
409: for (c = cStart; c < cEnd; ++c) {
410: const PetscInt *cone, *faceSizes, *faces;
411: const DMPolytopeType *faceTypes;
412: DMPolytopeType ct;
413: PetscInt numFaces, cf, foff = 0;
415: DMPlexGetCellType(dm, c, &ct);
416: DMPlexGetCone(dm, c, &cone);
417: DMPlexGetRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);
418: DMPlexSetCellType(idm, c, ct);
419: DMPlexSetConeSize(idm, c, numFaces);
420: for (cf = 0; cf < numFaces; foff += faceSizes[cf], ++cf) {
421: const PetscInt faceSize = faceSizes[cf];
422: const DMPolytopeType faceType = faceTypes[cf];
423: const PetscInt *face = &faces[foff];
424: PetscHashIJKLKey key;
425: PetscHashIter iter;
426: PetscBool missing;
427: PetscInt f;
429: if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not support faces of size %D > 4", faceSize);
430: key.i = face[0];
431: key.j = faceSize > 1 ? face[1] : PETSC_MAX_INT;
432: key.k = faceSize > 2 ? face[2] : PETSC_MAX_INT;
433: key.l = faceSize > 3 ? face[3] : PETSC_MAX_INT;
434: PetscSortInt(faceSize, (PetscInt *) &key);
435: PetscHashIJKLPut(faceTable, key, &iter, &missing);
436: if (missing) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Missing face (cell %D, lf %D)", c, cf);
437: PetscHashIJKLIterGet(faceTable, iter, &f);
438: DMPlexSetConeSize(idm, f, faceSize);
439: DMPlexSetCellType(idm, f, faceType);
440: }
441: DMPlexRestoreRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);
442: }
443: DMSetUp(idm);
444: /* Initialize cones so we do not need the bash table to tell us that a cone has been set */
445: {
446: PetscSection cs;
447: PetscInt *cones, csize;
449: DMPlexGetConeSection(idm, &cs);
450: DMPlexGetCones(idm, &cones);
451: PetscSectionGetStorageSize(cs, &csize);
452: for (c = 0; c < csize; ++c) cones[c] = -1;
453: }
454: /* Set cones */
455: for (d = 0; d <= depth; ++d) {
456: const PetscInt *cone;
457: PetscInt pStart, pEnd, p;
459: if (d == cellDepth) continue;
460: DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
461: for (p = pStart; p < pEnd; ++p) {
462: DMPlexGetCone(dm, p, &cone);
463: DMPlexSetCone(idm, p, cone);
464: DMPlexGetConeOrientation(dm, p, &cone);
465: DMPlexSetConeOrientation(idm, p, cone);
466: }
467: }
468: for (c = cStart; c < cEnd; ++c) {
469: const PetscInt *cone, *faceSizes, *faces;
470: const DMPolytopeType *faceTypes;
471: DMPolytopeType ct;
472: PetscInt numFaces, cf, foff = 0;
474: DMPlexGetCellType(dm, c, &ct);
475: DMPlexGetCone(dm, c, &cone);
476: DMPlexGetRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);
477: for (cf = 0; cf < numFaces; foff += faceSizes[cf], ++cf) {
478: DMPolytopeType faceType = faceTypes[cf];
479: const PetscInt faceSize = faceSizes[cf];
480: const PetscInt *face = &faces[foff];
481: const PetscInt *fcone;
482: PetscHashIJKLKey key;
483: PetscHashIter iter;
484: PetscBool missing;
485: PetscInt f;
487: if (faceSize > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Do not support faces of size %D > 4", faceSize);
488: key.i = face[0];
489: key.j = faceSize > 1 ? face[1] : PETSC_MAX_INT;
490: key.k = faceSize > 2 ? face[2] : PETSC_MAX_INT;
491: key.l = faceSize > 3 ? face[3] : PETSC_MAX_INT;
492: PetscSortInt(faceSize, (PetscInt *) &key);
493: PetscHashIJKLPut(faceTable, key, &iter, &missing);
494: PetscHashIJKLIterGet(faceTable, iter, &f);
495: DMPlexInsertCone(idm, c, cf, f);
496: DMPlexGetCone(idm, f, &fcone);
497: if (fcone[0] < 0) {DMPlexSetCone(idm, f, face);}
498: /* TODO This should be unnecessary since we have autoamtic orientation */
499: {
500: /* when matching hybrid faces in 3D, only few cases are possible.
501: Face traversal however can no longer follow the usual convention, this seems a serious issue to me */
502: PetscInt tquad_map[4][4] = { {PETSC_MIN_INT, 0,PETSC_MIN_INT,PETSC_MIN_INT},
503: { -1,PETSC_MIN_INT,PETSC_MIN_INT,PETSC_MIN_INT},
504: {PETSC_MIN_INT,PETSC_MIN_INT,PETSC_MIN_INT, 1},
505: {PETSC_MIN_INT,PETSC_MIN_INT, -2,PETSC_MIN_INT} };
506: PetscInt i, i2, j;
507: const PetscInt *cone;
508: PetscInt coneSize, ornt;
510: /* Orient face: Do not allow reverse orientation at the first vertex */
511: DMPlexGetConeSize(idm, f, &coneSize);
512: DMPlexGetCone(idm, f, &cone);
513: if (coneSize != faceSize) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of face vertices %D for face %D should be %D", coneSize, f, faceSize);
514: /* - First find the initial vertex */
515: for (i = 0; i < faceSize; ++i) if (face[0] == cone[i]) break;
516: /* If we want to compare tensor faces to regular faces, we have to flip them and take the else branch here */
517: if (faceType == DM_POLYTOPE_SEG_PRISM_TENSOR) {
518: /* find the second vertex */
519: for (i2 = 0; i2 < faceSize; ++i2) if (face[1] == cone[i2]) break;
520: switch (faceSize) {
521: case 2:
522: ornt = i ? -2 : 0;
523: break;
524: case 4:
525: ornt = tquad_map[i][i2];
526: break;
527: default: SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unhandled face size %D for face %D in cell %D", faceSize, f, c);
528: }
529: } else {
530: /* Try forward comparison */
531: for (j = 0; j < faceSize; ++j) if (face[j] != cone[(i+j)%faceSize]) break;
532: if (j == faceSize) {
533: if ((faceSize == 2) && (i == 1)) ornt = -2;
534: else ornt = i;
535: } else {
536: /* Try backward comparison */
537: for (j = 0; j < faceSize; ++j) if (face[j] != cone[(i+faceSize-j)%faceSize]) break;
538: if (j == faceSize) {
539: if (i == 0) ornt = -faceSize;
540: else ornt = -i;
541: } else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not determine orientation of face %D in cell %D", f, c);
542: }
543: }
544: DMPlexInsertConeOrientation(idm, c, cf, ornt);
545: }
546: }
547: DMPlexRestoreRawFaces_Internal(dm, ct, cone, &numFaces, &faceTypes, &faceSizes, &faces);
548: }
549: PetscHashIJKLDestroy(&faceTable);
550: DMPlexSymmetrize(idm);
551: DMPlexStratify(idm);
552: return(0);
553: }
555: static PetscErrorCode SortRmineRremoteByRemote_Private(PetscSF sf, PetscInt *rmine1[], PetscInt *rremote1[])
556: {
557: PetscInt nleaves;
558: PetscInt nranks;
559: const PetscMPIInt *ranks=NULL;
560: const PetscInt *roffset=NULL, *rmine=NULL, *rremote=NULL;
561: PetscInt n, o, r;
562: PetscErrorCode ierr;
565: PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, &rmine, &rremote);
566: nleaves = roffset[nranks];
567: PetscMalloc2(nleaves, rmine1, nleaves, rremote1);
568: for (r=0; r<nranks; r++) {
569: /* simultaneously sort rank-wise portions of rmine & rremote by values in rremote
570: - to unify order with the other side */
571: o = roffset[r];
572: n = roffset[r+1] - o;
573: PetscArraycpy(&(*rmine1)[o], &rmine[o], n);
574: PetscArraycpy(&(*rremote1)[o], &rremote[o], n);
575: PetscSortIntWithArray(n, &(*rremote1)[o], &(*rmine1)[o]);
576: }
577: return(0);
578: }
580: PetscErrorCode DMPlexOrientInterface_Internal(DM dm)
581: {
582: /* Here we only compare first 2 points of the cone. Full cone size would lead to stronger self-checking. */
583: PetscInt masterCone[2];
584: PetscInt (*roots)[2], (*leaves)[2];
585: PetscMPIInt (*rootsRanks)[2], (*leavesRanks)[2];
587: PetscSF sf=NULL;
588: const PetscInt *locals=NULL;
589: const PetscSFNode *remotes=NULL;
590: PetscInt nroots, nleaves, p, c;
591: PetscInt nranks, n, o, r;
592: const PetscMPIInt *ranks=NULL;
593: const PetscInt *roffset=NULL;
594: PetscInt *rmine1=NULL, *rremote1=NULL; /* rmine and rremote copies simultaneously sorted by rank and rremote */
595: const PetscInt *cone=NULL;
596: PetscInt coneSize, ind0;
597: MPI_Comm comm;
598: PetscMPIInt rank,size;
599: PetscInt debug = 0;
600: PetscErrorCode ierr;
603: DMGetPointSF(dm, &sf);
604: PetscSFGetGraph(sf, &nroots, &nleaves, &locals, &remotes);
605: if (nroots < 0) return(0);
606: PetscSFSetUp(sf);
607: PetscSFGetRootRanks(sf, &nranks, &ranks, &roffset, NULL, NULL);
608: DMViewFromOptions(dm, NULL, "-before_fix_dm_view");
609: if (PetscDefined(USE_DEBUG)) {DMPlexCheckPointSF(dm);}
610: SortRmineRremoteByRemote_Private(sf, &rmine1, &rremote1);
611: PetscMalloc4(nroots, &roots, nroots, &leaves, nroots, &rootsRanks, nroots, &leavesRanks);
612: PetscObjectGetComm((PetscObject) dm, &comm);
613: MPI_Comm_rank(comm, &rank);
614: MPI_Comm_size(comm, &size);
615: for (p = 0; p < nroots; ++p) {
616: DMPlexGetConeSize(dm, p, &coneSize);
617: DMPlexGetCone(dm, p, &cone);
618: if (coneSize < 2) {
619: for (c = 0; c < 2; c++) {
620: roots[p][c] = -1;
621: rootsRanks[p][c] = -1;
622: }
623: continue;
624: }
625: /* Translate all points to root numbering */
626: for (c = 0; c < 2; c++) {
627: PetscFindInt(cone[c], nleaves, locals, &ind0);
628: if (ind0 < 0) {
629: roots[p][c] = cone[c];
630: rootsRanks[p][c] = rank;
631: } else {
632: roots[p][c] = remotes[ind0].index;
633: rootsRanks[p][c] = remotes[ind0].rank;
634: }
635: }
636: }
637: for (p = 0; p < nroots; ++p) {
638: for (c = 0; c < 2; c++) {
639: leaves[p][c] = -2;
640: leavesRanks[p][c] = -2;
641: }
642: }
643: PetscSFBcastBegin(sf, MPIU_2INT, roots, leaves);
644: PetscSFBcastBegin(sf, MPI_2INT, rootsRanks, leavesRanks);
645: PetscSFBcastEnd(sf, MPIU_2INT, roots, leaves);
646: PetscSFBcastEnd(sf, MPI_2INT, rootsRanks, leavesRanks);
647: if (debug) {PetscSynchronizedFlush(comm, NULL);}
648: if (debug && rank == 0) {PetscSynchronizedPrintf(comm, "Referenced roots\n");}
649: for (p = 0; p < nroots; ++p) {
650: if (leaves[p][0] < 0) continue;
651: DMPlexGetConeSize(dm, p, &coneSize);
652: DMPlexGetCone(dm, p, &cone);
653: if (debug) {PetscSynchronizedPrintf(comm, "[%d] %4D: cone=[%4D %4D] roots=[(%d,%4D) (%d,%4D)] leaves=[(%d,%4D) (%d,%4D)]", rank, p, cone[0], cone[1], rootsRanks[p][0], roots[p][0], rootsRanks[p][1], roots[p][1], leavesRanks[p][0], leaves[p][0], leavesRanks[p][1], leaves[p][1]);}
654: if ((leaves[p][0] != roots[p][0]) || (leaves[p][1] != roots[p][1]) || (leavesRanks[p][0] != rootsRanks[p][0]) || (leavesRanks[p][1] != rootsRanks[p][1])) {
655: /* Translate these two leaves to my cone points; masterCone means desired order p's cone points */
656: for (c = 0; c < 2; c++) {
657: if (leavesRanks[p][c] == rank) {
658: /* A local leave is just taken as it is */
659: masterCone[c] = leaves[p][c];
660: continue;
661: }
662: /* Find index of rank leavesRanks[p][c] among remote ranks */
663: /* No need for PetscMPIIntCast because these integers were originally cast from PetscMPIInt. */
664: PetscFindMPIInt((PetscMPIInt)leavesRanks[p][c], nranks, ranks, &r);
665: if (PetscUnlikely(r < 0)) SETERRQ7(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D cone[%D]=%D root (%d,%D) leave (%d,%D): leave rank not found among remote ranks",p,c,cone[c],rootsRanks[p][c],roots[p][c],leavesRanks[p][c],leaves[p][c]);
666: if (PetscUnlikely(ranks[r] < 0 || ranks[r] >= size)) SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_PLIB, "p=%D c=%D commsize=%d: ranks[%D] = %d makes no sense",p,c,size,r,ranks[r]);
667: /* Find point leaves[p][c] among remote points aimed at rank leavesRanks[p][c] */
668: o = roffset[r];
669: n = roffset[r+1] - o;
670: PetscFindInt(leaves[p][c], n, &rremote1[o], &ind0);
671: if (PetscUnlikely(ind0 < 0)) SETERRQ7(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Point %D cone[%D]=%D root (%d,%D) leave (%d,%D): corresponding remote point not found - it seems there is missing connection in point SF!",p,c,cone[c],rootsRanks[p][c],roots[p][c],leavesRanks[p][c],leaves[p][c]);
672: /* Get the corresponding local point */
673: masterCone[c] = rmine1[o+ind0];
674: }
675: if (debug) {PetscSynchronizedPrintf(comm, " masterCone=[%4D %4D]\n", masterCone[0], masterCone[1]);}
676: /* Set the desired order of p's cone points and fix orientations accordingly */
677: /* Vaclav's note: Here we only compare first 2 points of the cone. Full cone size would lead to stronger self-checking. */
678: DMPlexOrientCell(dm, p, 2, masterCone);
679: } else if (debug) {PetscSynchronizedPrintf(comm, " ==\n");}
680: }
681: if (debug) {
682: PetscSynchronizedFlush(comm, NULL);
683: MPI_Barrier(comm);
684: }
685: DMViewFromOptions(dm, NULL, "-after_fix_dm_view");
686: PetscFree4(roots, leaves, rootsRanks, leavesRanks);
687: PetscFree2(rmine1, rremote1);
688: return(0);
689: }
691: static PetscErrorCode IntArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], const char valname[], PetscInt n, const PetscInt a[])
692: {
693: PetscInt idx;
694: PetscMPIInt rank;
695: PetscBool flg;
699: PetscOptionsHasName(NULL, NULL, opt, &flg);
700: if (!flg) return(0);
701: MPI_Comm_rank(comm, &rank);
702: PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);
703: for (idx = 0; idx < n; ++idx) {PetscSynchronizedPrintf(comm, "[%d]%s %D %s %D\n", rank, idxname, idx, valname, a[idx]);}
704: PetscSynchronizedFlush(comm, NULL);
705: return(0);
706: }
708: static PetscErrorCode SFNodeArrayViewFromOptions(MPI_Comm comm, const char opt[], const char name[], const char idxname[], PetscInt n, const PetscSFNode a[])
709: {
710: PetscInt idx;
711: PetscMPIInt rank;
712: PetscBool flg;
716: PetscOptionsHasName(NULL, NULL, opt, &flg);
717: if (!flg) return(0);
718: MPI_Comm_rank(comm, &rank);
719: PetscSynchronizedPrintf(comm, "[%d]%s:\n", rank, name);
720: if (idxname) {
721: for (idx = 0; idx < n; ++idx) {PetscSynchronizedPrintf(comm, "[%d]%s %D rank %D index %D\n", rank, idxname, idx, a[idx].rank, a[idx].index);}
722: } else {
723: for (idx = 0; idx < n; ++idx) {PetscSynchronizedPrintf(comm, "[%d]rank %D index %D\n", rank, a[idx].rank, a[idx].index);}
724: }
725: PetscSynchronizedFlush(comm, NULL);
726: return(0);
727: }
729: static PetscErrorCode DMPlexMapToLocalPoint(DM dm, PetscHMapIJ remotehash, PetscSFNode remotePoint, PetscInt *localPoint)
730: {
731: PetscSF sf;
732: const PetscInt *locals;
733: PetscMPIInt rank;
734: PetscErrorCode ierr;
737: MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
738: DMGetPointSF(dm, &sf);
739: PetscSFGetGraph(sf, NULL, NULL, &locals, NULL);
740: if (remotePoint.rank == rank) {
741: *localPoint = remotePoint.index;
742: } else {
743: PetscHashIJKey key;
744: PetscInt l;
746: key.i = remotePoint.index;
747: key.j = remotePoint.rank;
748: PetscHMapIJGet(remotehash, key, &l);
749: if (l >= 0) {
750: *localPoint = locals[l];
751: } else PetscFunctionReturn(1);
752: }
753: return(0);
754: }
756: static PetscErrorCode DMPlexMapToGlobalPoint(DM dm, PetscInt localPoint, PetscSFNode *remotePoint)
757: {
758: PetscSF sf;
759: const PetscInt *locals, *rootdegree;
760: const PetscSFNode *remotes;
761: PetscInt Nl, l;
762: PetscMPIInt rank;
763: PetscErrorCode ierr;
766: MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
767: DMGetPointSF(dm, &sf);
768: PetscSFGetGraph(sf, NULL, &Nl, &locals, &remotes);
769: if (Nl < 0) goto owned;
770: PetscSFComputeDegreeBegin(sf, &rootdegree);
771: PetscSFComputeDegreeEnd(sf, &rootdegree);
772: if (rootdegree[localPoint]) goto owned;
773: PetscFindInt(localPoint, Nl, locals, &l);
774: if (l < 0) PetscFunctionReturn(1);
775: *remotePoint = remotes[l];
776: return(0);
777: owned:
778: remotePoint->rank = rank;
779: remotePoint->index = localPoint;
780: return(0);
781: }
784: static PetscErrorCode DMPlexPointIsShared(DM dm, PetscInt p, PetscBool *isShared)
785: {
786: PetscSF sf;
787: const PetscInt *locals, *rootdegree;
788: PetscInt Nl, idx;
789: PetscErrorCode ierr;
792: *isShared = PETSC_FALSE;
793: DMGetPointSF(dm, &sf);
794: PetscSFGetGraph(sf, NULL, &Nl, &locals, NULL);
795: if (Nl < 0) return(0);
796: PetscFindInt(p, Nl, locals, &idx);
797: if (idx >= 0) {*isShared = PETSC_TRUE; return(0);}
798: PetscSFComputeDegreeBegin(sf, &rootdegree);
799: PetscSFComputeDegreeEnd(sf, &rootdegree);
800: if (rootdegree[p] > 0) *isShared = PETSC_TRUE;
801: return(0);
802: }
804: static PetscErrorCode DMPlexConeIsShared(DM dm, PetscInt p, PetscBool *isShared)
805: {
806: const PetscInt *cone;
807: PetscInt coneSize, c;
808: PetscBool cShared = PETSC_TRUE;
809: PetscErrorCode ierr;
812: DMPlexGetConeSize(dm, p, &coneSize);
813: DMPlexGetCone(dm, p, &cone);
814: for (c = 0; c < coneSize; ++c) {
815: PetscBool pointShared;
817: DMPlexPointIsShared(dm, cone[c], &pointShared);
818: cShared = (PetscBool) (cShared && pointShared);
819: }
820: *isShared = coneSize ? cShared : PETSC_FALSE;
821: return(0);
822: }
824: static PetscErrorCode DMPlexGetConeMinimum(DM dm, PetscInt p, PetscSFNode *cpmin)
825: {
826: const PetscInt *cone;
827: PetscInt coneSize, c;
828: PetscSFNode cmin = {PETSC_MAX_INT, PETSC_MAX_INT}, missing = {-1, -1};
829: PetscErrorCode ierr;
832: DMPlexGetConeSize(dm, p, &coneSize);
833: DMPlexGetCone(dm, p, &cone);
834: for (c = 0; c < coneSize; ++c) {
835: PetscSFNode rcp;
837: DMPlexMapToGlobalPoint(dm, cone[c], &rcp);
838: if (ierr) {
839: cmin = missing;
840: } else {
841: cmin = (rcp.rank < cmin.rank) || (rcp.rank == cmin.rank && rcp.index < cmin.index) ? rcp : cmin;
842: }
843: }
844: *cpmin = coneSize ? cmin : missing;
845: return(0);
846: }
848: /*
849: Each shared face has an entry in the candidates array:
850: (-1, coneSize-1), {(global cone point)}
851: where the set is missing the point p which we use as the key for the face
852: */
853: static PetscErrorCode DMPlexAddSharedFace_Private(DM dm, PetscSection candidateSection, PetscSFNode candidates[], PetscHMapIJ faceHash, PetscInt p, PetscBool debug)
854: {
855: MPI_Comm comm;
856: const PetscInt *support;
857: PetscInt supportSize, s, off = 0, idx = 0, overlap, cellHeight, height;
858: PetscMPIInt rank;
859: PetscErrorCode ierr;
862: PetscObjectGetComm((PetscObject) dm, &comm);
863: MPI_Comm_rank(comm, &rank);
864: DMPlexGetOverlap(dm, &overlap);
865: DMPlexGetVTKCellHeight(dm, &cellHeight);
866: DMPlexGetPointHeight(dm, p, &height);
867: if (!overlap && height <= cellHeight+1) {
868: /* cells can't be shared for non-overlapping meshes */
869: if (debug) {PetscSynchronizedPrintf(comm, "[%d] Skipping face %D to avoid adding cell to hashmap since this is nonoverlapping mesh\n", rank, p);}
870: return(0);
871: }
872: DMPlexGetSupportSize(dm, p, &supportSize);
873: DMPlexGetSupport(dm, p, &support);
874: if (candidates) {PetscSectionGetOffset(candidateSection, p, &off);}
875: for (s = 0; s < supportSize; ++s) {
876: const PetscInt face = support[s];
877: const PetscInt *cone;
878: PetscSFNode cpmin={-1,-1}, rp={-1,-1};
879: PetscInt coneSize, c, f;
880: PetscBool isShared = PETSC_FALSE;
881: PetscHashIJKey key;
883: /* Only add point once */
884: if (debug) {PetscSynchronizedPrintf(comm, "[%d] Support face %D\n", rank, face);}
885: key.i = p;
886: key.j = face;
887: PetscHMapIJGet(faceHash, key, &f);
888: if (f >= 0) continue;
889: DMPlexConeIsShared(dm, face, &isShared);
890: DMPlexGetConeMinimum(dm, face, &cpmin);
891: DMPlexMapToGlobalPoint(dm, p, &rp);
892: if (debug) {
893: PetscSynchronizedPrintf(comm, "[%d] Face point %D is shared: %d\n", rank, face, (int) isShared);
894: PetscSynchronizedPrintf(comm, "[%d] Global point (%D, %D) Min Cone Point (%D, %D)\n", rank, rp.rank, rp.index, cpmin.rank, cpmin.index);
895: }
896: if (isShared && (rp.rank == cpmin.rank && rp.index == cpmin.index)) {
897: PetscHMapIJSet(faceHash, key, p);
898: if (candidates) {
899: if (debug) {PetscSynchronizedPrintf(comm, "[%d] Adding shared face %D at idx %D\n[%d] ", rank, face, idx, rank);}
900: DMPlexGetConeSize(dm, face, &coneSize);
901: DMPlexGetCone(dm, face, &cone);
902: candidates[off+idx].rank = -1;
903: candidates[off+idx++].index = coneSize-1;
904: candidates[off+idx].rank = rank;
905: candidates[off+idx++].index = face;
906: for (c = 0; c < coneSize; ++c) {
907: const PetscInt cp = cone[c];
909: if (cp == p) continue;
910: DMPlexMapToGlobalPoint(dm, cp, &candidates[off+idx]);
911: if (debug) {PetscSynchronizedPrintf(comm, " (%D,%D)", candidates[off+idx].rank, candidates[off+idx].index);}
912: ++idx;
913: }
914: if (debug) {PetscSynchronizedPrintf(comm, "\n");}
915: } else {
916: /* Add cone size to section */
917: if (debug) {PetscSynchronizedPrintf(comm, "[%d] Scheduling shared face %D\n", rank, face);}
918: DMPlexGetConeSize(dm, face, &coneSize);
919: PetscHMapIJSet(faceHash, key, p);
920: PetscSectionAddDof(candidateSection, p, coneSize+1);
921: }
922: }
923: }
924: return(0);
925: }
927: /*@
928: DMPlexInterpolatePointSF - Insert interpolated points in the overlap into the PointSF in parallel, following local interpolation
930: Collective on dm
932: Input Parameters:
933: + dm - The interpolated DM
934: - pointSF - The initial SF without interpolated points
936: Output Parameter:
937: . pointSF - The SF including interpolated points
939: Level: developer
941: Note: All debugging for this process can be turned on with the options: -dm_interp_pre_view -petscsf_interp_pre_view -petscsection_interp_candidate_view -petscsection_interp_candidate_remote_view -petscsection_interp_claim_view -petscsf_interp_pre_view -dmplex_interp_debug
943: .seealso: DMPlexInterpolate(), DMPlexUninterpolate()
944: @*/
945: PetscErrorCode DMPlexInterpolatePointSF(DM dm, PetscSF pointSF)
946: {
947: MPI_Comm comm;
948: PetscHMapIJ remoteHash;
949: PetscHMapI claimshash;
950: PetscSection candidateSection, candidateRemoteSection, claimSection;
951: PetscSFNode *candidates, *candidatesRemote, *claims;
952: const PetscInt *localPoints, *rootdegree;
953: const PetscSFNode *remotePoints;
954: PetscInt ov, Nr, r, Nl, l;
955: PetscInt candidatesSize, candidatesRemoteSize, claimsSize;
956: PetscBool flg, debug = PETSC_FALSE;
957: PetscMPIInt rank;
958: PetscErrorCode ierr;
963: DMPlexIsDistributed(dm, &flg);
964: if (!flg) return(0);
965: /* Set initial SF so that lower level queries work */
966: DMSetPointSF(dm, pointSF);
967: PetscObjectGetComm((PetscObject) dm, &comm);
968: MPI_Comm_rank(comm, &rank);
969: DMPlexGetOverlap(dm, &ov);
970: if (ov) SETERRQ(comm, PETSC_ERR_SUP, "Interpolation of overlapped DMPlex not implemented yet");
971: PetscOptionsHasName(NULL, ((PetscObject) dm)->prefix, "-dmplex_interp_debug", &debug);
972: PetscObjectViewFromOptions((PetscObject) dm, NULL, "-dm_interp_pre_view");
973: PetscObjectViewFromOptions((PetscObject) pointSF, NULL, "-petscsf_interp_pre_view");
974: PetscLogEventBegin(DMPLEX_InterpolateSF,dm,0,0,0);
975: /* Step 0: Precalculations */
976: PetscSFGetGraph(pointSF, &Nr, &Nl, &localPoints, &remotePoints);
977: if (Nr < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "This DMPlex is distributed but input PointSF has no graph set");
978: PetscHMapIJCreate(&remoteHash);
979: for (l = 0; l < Nl; ++l) {
980: PetscHashIJKey key;
981: key.i = remotePoints[l].index;
982: key.j = remotePoints[l].rank;
983: PetscHMapIJSet(remoteHash, key, l);
984: }
985: /* Compute root degree to identify shared points */
986: PetscSFComputeDegreeBegin(pointSF, &rootdegree);
987: PetscSFComputeDegreeEnd(pointSF, &rootdegree);
988: IntArrayViewFromOptions(comm, "-interp_root_degree_view", "Root degree", "point", "degree", Nr, rootdegree);
989: /*
990: 1) Loop over each leaf point $p$ at depth $d$ in the SF
991: \item Get set $F(p)$ of faces $f$ in the support of $p$ for which
992: \begin{itemize}
993: \item all cone points of $f$ are shared
994: \item $p$ is the cone point with smallest canonical number
995: \end{itemize}
996: \item Send $F(p)$ and the cone of each face to the active root point $r(p)$
997: \item At the root, if at least two faces with a given cone are present, including a local face, mark the face as shared \label{alg:rootStep} and choose the root face
998: \item Send the root face from the root back to all leaf process
999: \item Leaf processes add the shared face to the SF
1000: */
1001: /* Step 1: Construct section+SFNode array
1002: The section has entries for all shared faces for which we have a leaf point in the cone
1003: The array holds candidate shared faces, each face is refered to by the leaf point */
1004: PetscSectionCreate(comm, &candidateSection);
1005: PetscSectionSetChart(candidateSection, 0, Nr);
1006: {
1007: PetscHMapIJ faceHash;
1009: PetscHMapIJCreate(&faceHash);
1010: for (l = 0; l < Nl; ++l) {
1011: const PetscInt p = localPoints[l];
1013: if (debug) {PetscSynchronizedPrintf(comm, "[%d] First pass leaf point %D\n", rank, p);}
1014: DMPlexAddSharedFace_Private(dm, candidateSection, NULL, faceHash, p, debug);
1015: }
1016: PetscHMapIJClear(faceHash);
1017: PetscSectionSetUp(candidateSection);
1018: PetscSectionGetStorageSize(candidateSection, &candidatesSize);
1019: PetscMalloc1(candidatesSize, &candidates);
1020: for (l = 0; l < Nl; ++l) {
1021: const PetscInt p = localPoints[l];
1023: if (debug) {PetscSynchronizedPrintf(comm, "[%d] Second pass leaf point %D\n", rank, p);}
1024: DMPlexAddSharedFace_Private(dm, candidateSection, candidates, faceHash, p, debug);
1025: }
1026: PetscHMapIJDestroy(&faceHash);
1027: if (debug) {PetscSynchronizedFlush(comm, NULL);}
1028: }
1029: PetscObjectSetName((PetscObject) candidateSection, "Candidate Section");
1030: PetscObjectViewFromOptions((PetscObject) candidateSection, NULL, "-petscsection_interp_candidate_view");
1031: SFNodeArrayViewFromOptions(comm, "-petscsection_interp_candidate_view", "Candidates", NULL, candidatesSize, candidates);
1032: /* Step 2: Gather candidate section / array pair into the root partition via inverse(multi(pointSF)). */
1033: /* Note that this section is indexed by offsets into leaves, not by point number */
1034: {
1035: PetscSF sfMulti, sfInverse, sfCandidates;
1036: PetscInt *remoteOffsets;
1038: PetscSFGetMultiSF(pointSF, &sfMulti);
1039: PetscSFCreateInverseSF(sfMulti, &sfInverse);
1040: PetscSectionCreate(comm, &candidateRemoteSection);
1041: PetscSFDistributeSection(sfInverse, candidateSection, &remoteOffsets, candidateRemoteSection);
1042: PetscSFCreateSectionSF(sfInverse, candidateSection, remoteOffsets, candidateRemoteSection, &sfCandidates);
1043: PetscSectionGetStorageSize(candidateRemoteSection, &candidatesRemoteSize);
1044: PetscMalloc1(candidatesRemoteSize, &candidatesRemote);
1045: PetscSFBcastBegin(sfCandidates, MPIU_2INT, candidates, candidatesRemote);
1046: PetscSFBcastEnd(sfCandidates, MPIU_2INT, candidates, candidatesRemote);
1047: PetscSFDestroy(&sfInverse);
1048: PetscSFDestroy(&sfCandidates);
1049: PetscFree(remoteOffsets);
1051: PetscObjectSetName((PetscObject) candidateRemoteSection, "Remote Candidate Section");
1052: PetscObjectViewFromOptions((PetscObject) candidateRemoteSection, NULL, "-petscsection_interp_candidate_remote_view");
1053: SFNodeArrayViewFromOptions(comm, "-petscsection_interp_candidate_remote_view", "Remote Candidates", NULL, candidatesRemoteSize, candidatesRemote);
1054: }
1055: /* Step 3: At the root, if at least two faces with a given cone are present, including a local face, mark the face as shared and choose the root face */
1056: {
1057: PetscHashIJKLRemote faceTable;
1058: PetscInt idx, idx2;
1060: PetscHashIJKLRemoteCreate(&faceTable);
1061: /* There is a section point for every leaf attached to a given root point */
1062: for (r = 0, idx = 0, idx2 = 0; r < Nr; ++r) {
1063: PetscInt deg;
1065: for (deg = 0; deg < rootdegree[r]; ++deg, ++idx) {
1066: PetscInt offset, dof, d;
1068: PetscSectionGetDof(candidateRemoteSection, idx, &dof);
1069: PetscSectionGetOffset(candidateRemoteSection, idx, &offset);
1070: /* dof may include many faces from the remote process */
1071: for (d = 0; d < dof; ++d) {
1072: const PetscInt hidx = offset+d;
1073: const PetscInt Np = candidatesRemote[hidx].index+1;
1074: const PetscSFNode rface = candidatesRemote[hidx+1];
1075: const PetscSFNode *fcone = &candidatesRemote[hidx+2];
1076: PetscSFNode fcp0;
1077: const PetscSFNode pmax = {PETSC_MAX_INT, PETSC_MAX_INT};
1078: const PetscInt *join = NULL;
1079: PetscHashIJKLRemoteKey key;
1080: PetscHashIter iter;
1081: PetscBool missing;
1082: PetscInt points[1024], p, joinSize;
1084: if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking face (%D, %D) at (%D, %D, %D) with cone size %D\n", rank, rface.rank, rface.index, r, idx, d, Np);}
1085: if (Np > 4) SETERRQ6(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle face (%D, %D) at (%D, %D, %D) with %D cone points", rface.rank, rface.index, r, idx, d, Np);
1086: fcp0.rank = rank;
1087: fcp0.index = r;
1088: d += Np;
1089: /* Put remote face in hash table */
1090: key.i = fcp0;
1091: key.j = fcone[0];
1092: key.k = Np > 2 ? fcone[1] : pmax;
1093: key.l = Np > 3 ? fcone[2] : pmax;
1094: PetscSortSFNode(Np, (PetscSFNode *) &key);
1095: PetscHashIJKLRemotePut(faceTable, key, &iter, &missing);
1096: if (missing) {
1097: if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Setting remote face (%D, %D)\n", rank, rface.index, rface.rank);}
1098: PetscHashIJKLRemoteIterSet(faceTable, iter, rface);
1099: } else {
1100: PetscSFNode oface;
1102: PetscHashIJKLRemoteIterGet(faceTable, iter, &oface);
1103: if ((rface.rank < oface.rank) || (rface.rank == oface.rank && rface.index < oface.index)) {
1104: if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Replacing with remote face (%D, %D)\n", rank, rface.index, rface.rank);}
1105: PetscHashIJKLRemoteIterSet(faceTable, iter, rface);
1106: }
1107: }
1108: /* Check for local face */
1109: points[0] = r;
1110: for (p = 1; p < Np; ++p) {
1111: DMPlexMapToLocalPoint(dm, remoteHash, fcone[p-1], &points[p]);
1112: if (ierr) break; /* We got a point not in our overlap */
1113: if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Checking local candidate %D\n", rank, points[p]);}
1114: }
1115: if (ierr) continue;
1116: DMPlexGetJoin(dm, Np, points, &joinSize, &join);
1117: if (joinSize == 1) {
1118: PetscSFNode lface;
1119: PetscSFNode oface;
1121: /* Always replace with local face */
1122: lface.rank = rank;
1123: lface.index = join[0];
1124: PetscHashIJKLRemoteIterGet(faceTable, iter, &oface);
1125: if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Replacing (%D, %D) with local face (%D, %D)\n", rank, oface.index, oface.rank, lface.index, lface.rank);}
1126: PetscHashIJKLRemoteIterSet(faceTable, iter, lface);
1127: }
1128: DMPlexRestoreJoin(dm, Np, points, &joinSize, &join);
1129: }
1130: }
1131: /* Put back faces for this root */
1132: for (deg = 0; deg < rootdegree[r]; ++deg, ++idx2) {
1133: PetscInt offset, dof, d;
1135: PetscSectionGetDof(candidateRemoteSection, idx2, &dof);
1136: PetscSectionGetOffset(candidateRemoteSection, idx2, &offset);
1137: /* dof may include many faces from the remote process */
1138: for (d = 0; d < dof; ++d) {
1139: const PetscInt hidx = offset+d;
1140: const PetscInt Np = candidatesRemote[hidx].index+1;
1141: const PetscSFNode *fcone = &candidatesRemote[hidx+2];
1142: PetscSFNode fcp0;
1143: const PetscSFNode pmax = {PETSC_MAX_INT, PETSC_MAX_INT};
1144: PetscHashIJKLRemoteKey key;
1145: PetscHashIter iter;
1146: PetscBool missing;
1148: if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Entering face at (%D, %D)\n", rank, r, idx);}
1149: if (Np > 4) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle faces with %D cone points", Np);
1150: fcp0.rank = rank;
1151: fcp0.index = r;
1152: d += Np;
1153: /* Find remote face in hash table */
1154: key.i = fcp0;
1155: key.j = fcone[0];
1156: key.k = Np > 2 ? fcone[1] : pmax;
1157: key.l = Np > 3 ? fcone[2] : pmax;
1158: PetscSortSFNode(Np, (PetscSFNode *) &key);
1159: if (debug) {PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] key (%D, %D) (%D, %D) (%D, %D) (%D, %D)\n", rank, key.i.rank, key.i.index, key.j.rank, key.j.index, key.k.rank, key.k.index, key.l.rank, key.l.index);}
1160: PetscHashIJKLRemotePut(faceTable, key, &iter, &missing);
1161: if (missing) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Root %D Idx %D ought to have an associated face", r, idx2);
1162: else {PetscHashIJKLRemoteIterGet(faceTable, iter, &candidatesRemote[hidx]);}
1163: }
1164: }
1165: }
1166: if (debug) {PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), NULL);}
1167: PetscHashIJKLRemoteDestroy(&faceTable);
1168: }
1169: /* Step 4: Push back owned faces */
1170: {
1171: PetscSF sfMulti, sfClaims, sfPointNew;
1172: PetscSFNode *remotePointsNew;
1173: PetscInt *remoteOffsets, *localPointsNew;
1174: PetscInt pStart, pEnd, r, NlNew, p;
1176: /* 4) Push claims back to receiver via the MultiSF and derive new pointSF mapping on receiver */
1177: PetscSFGetMultiSF(pointSF, &sfMulti);
1178: PetscSectionCreate(comm, &claimSection);
1179: PetscSFDistributeSection(sfMulti, candidateRemoteSection, &remoteOffsets, claimSection);
1180: PetscSFCreateSectionSF(sfMulti, candidateRemoteSection, remoteOffsets, claimSection, &sfClaims);
1181: PetscSectionGetStorageSize(claimSection, &claimsSize);
1182: PetscMalloc1(claimsSize, &claims);
1183: for (p = 0; p < claimsSize; ++p) claims[p].rank = -1;
1184: PetscSFBcastBegin(sfClaims, MPIU_2INT, candidatesRemote, claims);
1185: PetscSFBcastEnd(sfClaims, MPIU_2INT, candidatesRemote, claims);
1186: PetscSFDestroy(&sfClaims);
1187: PetscFree(remoteOffsets);
1188: PetscObjectSetName((PetscObject) claimSection, "Claim Section");
1189: PetscObjectViewFromOptions((PetscObject) claimSection, NULL, "-petscsection_interp_claim_view");
1190: SFNodeArrayViewFromOptions(comm, "-petscsection_interp_claim_view", "Claims", NULL, claimsSize, claims);
1191: /* Step 5) Walk the original section of local supports and add an SF entry for each updated item */
1192: /* TODO I should not have to do a join here since I already put the face and its cone in the candidate section */
1193: PetscHMapICreate(&claimshash);
1194: for (r = 0; r < Nr; ++r) {
1195: PetscInt dof, off, d;
1197: if (debug) {PetscSynchronizedPrintf(comm, "[%d] Checking root for claims %D\n", rank, r);}
1198: PetscSectionGetDof(candidateSection, r, &dof);
1199: PetscSectionGetOffset(candidateSection, r, &off);
1200: for (d = 0; d < dof;) {
1201: if (claims[off+d].rank >= 0) {
1202: const PetscInt faceInd = off+d;
1203: const PetscInt Np = candidates[off+d].index;
1204: const PetscInt *join = NULL;
1205: PetscInt joinSize, points[1024], c;
1207: if (debug) {PetscSynchronizedPrintf(comm, "[%d] Found claim for remote point (%D, %D)\n", rank, claims[faceInd].rank, claims[faceInd].index);}
1208: points[0] = r;
1209: if (debug) {PetscSynchronizedPrintf(comm, "[%d] point %D\n", rank, points[0]);}
1210: for (c = 0, d += 2; c < Np; ++c, ++d) {
1211: DMPlexMapToLocalPoint(dm, remoteHash, candidates[off+d], &points[c+1]);
1212: if (debug) {PetscSynchronizedPrintf(comm, "[%d] point %D\n", rank, points[c+1]);}
1213: }
1214: DMPlexGetJoin(dm, Np+1, points, &joinSize, &join);
1215: if (joinSize == 1) {
1216: if (claims[faceInd].rank == rank) {
1217: if (debug) {PetscSynchronizedPrintf(comm, "[%d] Ignoring local face %D for non-remote partner\n", rank, join[0]);}
1218: } else {
1219: if (debug) {PetscSynchronizedPrintf(comm, "[%d] Found local face %D\n", rank, join[0]);}
1220: PetscHMapISet(claimshash, join[0], faceInd);
1221: }
1222: } else {
1223: if (debug) {PetscSynchronizedPrintf(comm, "[%d] Failed to find face\n", rank);}
1224: }
1225: DMPlexRestoreJoin(dm, Np+1, points, &joinSize, &join);
1226: } else {
1227: if (debug) {PetscSynchronizedPrintf(comm, "[%d] No claim for point %D\n", rank, r);}
1228: d += claims[off+d].index+1;
1229: }
1230: }
1231: }
1232: if (debug) {PetscSynchronizedFlush(comm, NULL);}
1233: /* Step 6) Create new pointSF from hashed claims */
1234: PetscHMapIGetSize(claimshash, &NlNew);
1235: DMPlexGetChart(dm, &pStart, &pEnd);
1236: PetscMalloc1(Nl + NlNew, &localPointsNew);
1237: PetscMalloc1(Nl + NlNew, &remotePointsNew);
1238: for (l = 0; l < Nl; ++l) {
1239: localPointsNew[l] = localPoints[l];
1240: remotePointsNew[l].index = remotePoints[l].index;
1241: remotePointsNew[l].rank = remotePoints[l].rank;
1242: }
1243: p = Nl;
1244: PetscHMapIGetKeys(claimshash, &p, localPointsNew);
1245: /* We sort new points, and assume they are numbered after all existing points */
1246: PetscSortInt(NlNew, &localPointsNew[Nl]);
1247: for (p = Nl; p < Nl + NlNew; ++p) {
1248: PetscInt off;
1249: PetscHMapIGet(claimshash, localPointsNew[p], &off);
1250: if (claims[off].rank < 0 || claims[off].index < 0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid claim for local point %D, (%D, %D)", localPointsNew[p], claims[off].rank, claims[off].index);
1251: remotePointsNew[p] = claims[off];
1252: }
1253: PetscSFCreate(comm, &sfPointNew);
1254: PetscSFSetGraph(sfPointNew, pEnd-pStart, Nl+NlNew, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);
1255: PetscSFSetUp(sfPointNew);
1256: DMSetPointSF(dm, sfPointNew);
1257: PetscObjectViewFromOptions((PetscObject) sfPointNew, NULL, "-petscsf_interp_view");
1258: PetscSFDestroy(&sfPointNew);
1259: PetscHMapIDestroy(&claimshash);
1260: }
1261: PetscHMapIJDestroy(&remoteHash);
1262: PetscSectionDestroy(&candidateSection);
1263: PetscSectionDestroy(&candidateRemoteSection);
1264: PetscSectionDestroy(&claimSection);
1265: PetscFree(candidates);
1266: PetscFree(candidatesRemote);
1267: PetscFree(claims);
1268: PetscLogEventEnd(DMPLEX_InterpolateSF,dm,0,0,0);
1269: return(0);
1270: }
1272: /*@
1273: DMPlexInterpolate - Take in a cell-vertex mesh and return one with all intermediate faces, edges, etc.
1275: Collective on dm
1277: Input Parameters:
1278: + dm - The DMPlex object with only cells and vertices
1279: - dmInt - The interpolated DM
1281: Output Parameter:
1282: . dmInt - The complete DMPlex object
1284: Level: intermediate
1286: Notes:
1287: It does not copy over the coordinates.
1289: Developer Notes:
1290: It sets plex->interpolated = DMPLEX_INTERPOLATED_FULL.
1292: .seealso: DMPlexUninterpolate(), DMPlexCreateFromCellListPetsc(), DMPlexCopyCoordinates()
1293: @*/
1294: PetscErrorCode DMPlexInterpolate(DM dm, DM *dmInt)
1295: {
1296: DMPlexInterpolatedFlag interpolated;
1297: DM idm, odm = dm;
1298: PetscSF sfPoint;
1299: PetscInt depth, dim, d;
1300: const char *name;
1301: PetscBool flg=PETSC_TRUE;
1307: PetscLogEventBegin(DMPLEX_Interpolate,dm,0,0,0);
1308: DMPlexGetDepth(dm, &depth);
1309: DMGetDimension(dm, &dim);
1310: DMPlexIsInterpolated(dm, &interpolated);
1311: if (interpolated == DMPLEX_INTERPOLATED_PARTIAL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not for partially interpolated meshes");
1312: if (interpolated == DMPLEX_INTERPOLATED_FULL) {
1313: PetscObjectReference((PetscObject) dm);
1314: idm = dm;
1315: } else {
1316: for (d = 1; d < dim; ++d) {
1317: /* Create interpolated mesh */
1318: DMCreate(PetscObjectComm((PetscObject)dm), &idm);
1319: DMSetType(idm, DMPLEX);
1320: DMSetDimension(idm, dim);
1321: if (depth > 0) {
1322: DMPlexInterpolateFaces_Internal(odm, 1, idm);
1323: DMGetPointSF(odm, &sfPoint);
1324: {
1325: /* TODO: We need to systematically fix cases of distributed Plexes with no graph set */
1326: PetscInt nroots;
1327: PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);
1328: if (nroots >= 0) {DMPlexInterpolatePointSF(idm, sfPoint);}
1329: }
1330: }
1331: if (odm != dm) {DMDestroy(&odm);}
1332: odm = idm;
1333: }
1334: PetscObjectGetName((PetscObject) dm, &name);
1335: PetscObjectSetName((PetscObject) idm, name);
1336: DMPlexCopyCoordinates(dm, idm);
1337: DMCopyLabels(dm, idm, PETSC_COPY_VALUES, PETSC_FALSE);
1338: PetscOptionsGetBool(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-dm_plex_interpolate_orient_interfaces", &flg, NULL);
1339: if (flg) {DMPlexOrientInterface_Internal(idm);}
1340: }
1341: {
1342: PetscBool isper;
1343: const PetscReal *maxCell, *L;
1344: const DMBoundaryType *bd;
1346: DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);
1347: DMSetPeriodicity(idm,isper,maxCell,L,bd);
1348: }
1349: /* This function makes the mesh fully interpolated on all ranks */
1350: {
1351: DM_Plex *plex = (DM_Plex *) idm->data;
1352: plex->interpolated = plex->interpolatedCollective = DMPLEX_INTERPOLATED_FULL;
1353: }
1354: *dmInt = idm;
1355: PetscLogEventEnd(DMPLEX_Interpolate,dm,0,0,0);
1356: return(0);
1357: }
1359: /*@
1360: DMPlexCopyCoordinates - Copy coordinates from one mesh to another with the same vertices
1362: Collective on dmA
1364: Input Parameter:
1365: . dmA - The DMPlex object with initial coordinates
1367: Output Parameter:
1368: . dmB - The DMPlex object with copied coordinates
1370: Level: intermediate
1372: Note: This is typically used when adding pieces other than vertices to a mesh
1374: .seealso: DMCopyLabels(), DMGetCoordinates(), DMGetCoordinatesLocal(), DMGetCoordinateDM(), DMGetCoordinateSection()
1375: @*/
1376: PetscErrorCode DMPlexCopyCoordinates(DM dmA, DM dmB)
1377: {
1378: Vec coordinatesA, coordinatesB;
1379: VecType vtype;
1380: PetscSection coordSectionA, coordSectionB;
1381: PetscScalar *coordsA, *coordsB;
1382: PetscInt spaceDim, Nf, vStartA, vStartB, vEndA, vEndB, coordSizeB, v, d;
1383: PetscInt cStartA, cEndA, cStartB, cEndB, cS, cE, cdim;
1384: PetscBool lc = PETSC_FALSE;
1390: if (dmA == dmB) return(0);
1391: DMGetCoordinateDim(dmA, &cdim);
1392: DMSetCoordinateDim(dmB, cdim);
1393: DMPlexGetDepthStratum(dmA, 0, &vStartA, &vEndA);
1394: DMPlexGetDepthStratum(dmB, 0, &vStartB, &vEndB);
1395: if ((vEndA-vStartA) != (vEndB-vStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of vertices in first DM %d != %d in the second DM", vEndA-vStartA, vEndB-vStartB);
1396: DMPlexGetHeightStratum(dmA, 0, &cStartA, &cEndA);
1397: DMPlexGetHeightStratum(dmB, 0, &cStartB, &cEndB);
1398: DMGetCoordinateSection(dmA, &coordSectionA);
1399: DMGetCoordinateSection(dmB, &coordSectionB);
1400: if (coordSectionA == coordSectionB) return(0);
1401: PetscSectionGetNumFields(coordSectionA, &Nf);
1402: if (!Nf) return(0);
1403: if (Nf > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of coordinate fields must be 1, not %D", Nf);
1404: if (!coordSectionB) {
1405: PetscInt dim;
1407: PetscSectionCreate(PetscObjectComm((PetscObject) coordSectionA), &coordSectionB);
1408: DMGetCoordinateDim(dmA, &dim);
1409: DMSetCoordinateSection(dmB, dim, coordSectionB);
1410: PetscObjectDereference((PetscObject) coordSectionB);
1411: }
1412: PetscSectionSetNumFields(coordSectionB, 1);
1413: PetscSectionGetFieldComponents(coordSectionA, 0, &spaceDim);
1414: PetscSectionSetFieldComponents(coordSectionB, 0, spaceDim);
1415: PetscSectionGetChart(coordSectionA, &cS, &cE);
1416: if (cStartA <= cS && cS < cEndA) { /* localized coordinates */
1417: if ((cEndA-cStartA) != (cEndB-cStartB)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "The number of cells in first DM %D != %D in the second DM", cEndA-cStartA, cEndB-cStartB);
1418: cS = cS - cStartA + cStartB;
1419: cE = vEndB;
1420: lc = PETSC_TRUE;
1421: } else {
1422: cS = vStartB;
1423: cE = vEndB;
1424: }
1425: PetscSectionSetChart(coordSectionB, cS, cE);
1426: for (v = vStartB; v < vEndB; ++v) {
1427: PetscSectionSetDof(coordSectionB, v, spaceDim);
1428: PetscSectionSetFieldDof(coordSectionB, v, 0, spaceDim);
1429: }
1430: if (lc) { /* localized coordinates */
1431: PetscInt c;
1433: for (c = cS-cStartB; c < cEndB-cStartB; c++) {
1434: PetscInt dof;
1436: PetscSectionGetDof(coordSectionA, c + cStartA, &dof);
1437: PetscSectionSetDof(coordSectionB, c + cStartB, dof);
1438: PetscSectionSetFieldDof(coordSectionB, c + cStartB, 0, dof);
1439: }
1440: }
1441: PetscSectionSetUp(coordSectionB);
1442: PetscSectionGetStorageSize(coordSectionB, &coordSizeB);
1443: DMGetCoordinatesLocal(dmA, &coordinatesA);
1444: VecCreate(PETSC_COMM_SELF, &coordinatesB);
1445: PetscObjectSetName((PetscObject) coordinatesB, "coordinates");
1446: VecSetSizes(coordinatesB, coordSizeB, PETSC_DETERMINE);
1447: VecGetBlockSize(coordinatesA, &d);
1448: VecSetBlockSize(coordinatesB, d);
1449: VecGetType(coordinatesA, &vtype);
1450: VecSetType(coordinatesB, vtype);
1451: VecGetArray(coordinatesA, &coordsA);
1452: VecGetArray(coordinatesB, &coordsB);
1453: for (v = 0; v < vEndB-vStartB; ++v) {
1454: PetscInt offA, offB;
1456: PetscSectionGetOffset(coordSectionA, v + vStartA, &offA);
1457: PetscSectionGetOffset(coordSectionB, v + vStartB, &offB);
1458: for (d = 0; d < spaceDim; ++d) {
1459: coordsB[offB+d] = coordsA[offA+d];
1460: }
1461: }
1462: if (lc) { /* localized coordinates */
1463: PetscInt c;
1465: for (c = cS-cStartB; c < cEndB-cStartB; c++) {
1466: PetscInt dof, offA, offB;
1468: PetscSectionGetOffset(coordSectionA, c + cStartA, &offA);
1469: PetscSectionGetOffset(coordSectionB, c + cStartB, &offB);
1470: PetscSectionGetDof(coordSectionA, c + cStartA, &dof);
1471: PetscArraycpy(coordsB + offB,coordsA + offA,dof);
1472: }
1473: }
1474: VecRestoreArray(coordinatesA, &coordsA);
1475: VecRestoreArray(coordinatesB, &coordsB);
1476: DMSetCoordinatesLocal(dmB, coordinatesB);
1477: VecDestroy(&coordinatesB);
1478: return(0);
1479: }
1481: /*@
1482: DMPlexUninterpolate - Take in a mesh with all intermediate faces, edges, etc. and return a cell-vertex mesh
1484: Collective on dm
1486: Input Parameter:
1487: . dm - The complete DMPlex object
1489: Output Parameter:
1490: . dmUnint - The DMPlex object with only cells and vertices
1492: Level: intermediate
1494: Notes:
1495: It does not copy over the coordinates.
1497: Developer Notes:
1498: It sets plex->interpolated = DMPLEX_INTERPOLATED_NONE.
1500: .seealso: DMPlexInterpolate(), DMPlexCreateFromCellListPetsc(), DMPlexCopyCoordinates()
1501: @*/
1502: PetscErrorCode DMPlexUninterpolate(DM dm, DM *dmUnint)
1503: {
1504: DMPlexInterpolatedFlag interpolated;
1505: DM udm;
1506: PetscInt dim, vStart, vEnd, cStart, cEnd, c, maxConeSize = 0, *cone;
1512: DMGetDimension(dm, &dim);
1513: DMPlexIsInterpolated(dm, &interpolated);
1514: if (interpolated == DMPLEX_INTERPOLATED_PARTIAL) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not for partially interpolated meshes");
1515: if (interpolated == DMPLEX_INTERPOLATED_NONE || dim <= 1) {
1516: /* in case dim <= 1 just keep the DMPLEX_INTERPOLATED_FULL flag */
1517: PetscObjectReference((PetscObject) dm);
1518: *dmUnint = dm;
1519: return(0);
1520: }
1521: DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);
1522: DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);
1523: DMCreate(PetscObjectComm((PetscObject) dm), &udm);
1524: DMSetType(udm, DMPLEX);
1525: DMSetDimension(udm, dim);
1526: DMPlexSetChart(udm, cStart, vEnd);
1527: for (c = cStart; c < cEnd; ++c) {
1528: PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
1530: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
1531: for (cl = 0; cl < closureSize*2; cl += 2) {
1532: const PetscInt p = closure[cl];
1534: if ((p >= vStart) && (p < vEnd)) ++coneSize;
1535: }
1536: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
1537: DMPlexSetConeSize(udm, c, coneSize);
1538: maxConeSize = PetscMax(maxConeSize, coneSize);
1539: }
1540: DMSetUp(udm);
1541: PetscMalloc1(maxConeSize, &cone);
1542: for (c = cStart; c < cEnd; ++c) {
1543: PetscInt *closure = NULL, closureSize, cl, coneSize = 0;
1545: DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
1546: for (cl = 0; cl < closureSize*2; cl += 2) {
1547: const PetscInt p = closure[cl];
1549: if ((p >= vStart) && (p < vEnd)) cone[coneSize++] = p;
1550: }
1551: DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure);
1552: DMPlexSetCone(udm, c, cone);
1553: }
1554: PetscFree(cone);
1555: DMPlexSymmetrize(udm);
1556: DMPlexStratify(udm);
1557: /* Reduce SF */
1558: {
1559: PetscSF sfPoint, sfPointUn;
1560: const PetscSFNode *remotePoints;
1561: const PetscInt *localPoints;
1562: PetscSFNode *remotePointsUn;
1563: PetscInt *localPointsUn;
1564: PetscInt vEnd, numRoots, numLeaves, l;
1565: PetscInt numLeavesUn = 0, n = 0;
1566: PetscErrorCode ierr;
1568: /* Get original SF information */
1569: DMGetPointSF(dm, &sfPoint);
1570: DMGetPointSF(udm, &sfPointUn);
1571: DMPlexGetDepthStratum(dm, 0, NULL, &vEnd);
1572: PetscSFGetGraph(sfPoint, &numRoots, &numLeaves, &localPoints, &remotePoints);
1573: /* Allocate space for cells and vertices */
1574: for (l = 0; l < numLeaves; ++l) if (localPoints[l] < vEnd) numLeavesUn++;
1575: /* Fill in leaves */
1576: if (vEnd >= 0) {
1577: PetscMalloc1(numLeavesUn, &remotePointsUn);
1578: PetscMalloc1(numLeavesUn, &localPointsUn);
1579: for (l = 0; l < numLeaves; l++) {
1580: if (localPoints[l] < vEnd) {
1581: localPointsUn[n] = localPoints[l];
1582: remotePointsUn[n].rank = remotePoints[l].rank;
1583: remotePointsUn[n].index = remotePoints[l].index;
1584: ++n;
1585: }
1586: }
1587: if (n != numLeavesUn) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistent number of leaves %d != %d", n, numLeavesUn);
1588: PetscSFSetGraph(sfPointUn, vEnd, numLeavesUn, localPointsUn, PETSC_OWN_POINTER, remotePointsUn, PETSC_OWN_POINTER);
1589: }
1590: }
1591: {
1592: PetscBool isper;
1593: const PetscReal *maxCell, *L;
1594: const DMBoundaryType *bd;
1596: DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);
1597: DMSetPeriodicity(udm,isper,maxCell,L,bd);
1598: }
1599: /* This function makes the mesh fully uninterpolated on all ranks */
1600: {
1601: DM_Plex *plex = (DM_Plex *) udm->data;
1602: plex->interpolated = plex->interpolatedCollective = DMPLEX_INTERPOLATED_NONE;
1603: }
1604: *dmUnint = udm;
1605: return(0);
1606: }
1608: static PetscErrorCode DMPlexIsInterpolated_Internal(DM dm, DMPlexInterpolatedFlag *interpolated)
1609: {
1610: PetscInt coneSize, depth, dim, h, p, pStart, pEnd;
1611: MPI_Comm comm;
1615: PetscObjectGetComm((PetscObject)dm, &comm);
1616: DMPlexGetDepth(dm, &depth);
1617: DMGetDimension(dm, &dim);
1619: if (depth == dim) {
1620: *interpolated = DMPLEX_INTERPOLATED_FULL;
1621: if (!dim) goto finish;
1623: /* Check points at height = dim are vertices (have no cones) */
1624: DMPlexGetHeightStratum(dm, dim, &pStart, &pEnd);
1625: for (p=pStart; p<pEnd; p++) {
1626: DMPlexGetConeSize(dm, p, &coneSize);
1627: if (coneSize) {
1628: *interpolated = DMPLEX_INTERPOLATED_PARTIAL;
1629: goto finish;
1630: }
1631: }
1633: /* Check points at height < dim have cones */
1634: for (h=0; h<dim; h++) {
1635: DMPlexGetHeightStratum(dm, h, &pStart, &pEnd);
1636: for (p=pStart; p<pEnd; p++) {
1637: DMPlexGetConeSize(dm, p, &coneSize);
1638: if (!coneSize) {
1639: *interpolated = DMPLEX_INTERPOLATED_PARTIAL;
1640: goto finish;
1641: }
1642: }
1643: }
1644: } else if (depth == 1) {
1645: *interpolated = DMPLEX_INTERPOLATED_NONE;
1646: } else {
1647: *interpolated = DMPLEX_INTERPOLATED_PARTIAL;
1648: }
1649: finish:
1650: return(0);
1651: }
1653: /*@
1654: DMPlexIsInterpolated - Find out to what extent the DMPlex is topologically interpolated.
1656: Not Collective
1658: Input Parameter:
1659: . dm - The DM object
1661: Output Parameter:
1662: . interpolated - Flag whether the DM is interpolated
1664: Level: intermediate
1666: Notes:
1667: Unlike DMPlexIsInterpolatedCollective(), this is NOT collective
1668: so the results can be different on different ranks in special cases.
1669: However, DMPlexInterpolate() guarantees the result is the same on all.
1671: Unlike DMPlexIsInterpolatedCollective(), this cannot return DMPLEX_INTERPOLATED_MIXED.
1673: Developer Notes:
1674: Initially, plex->interpolated = DMPLEX_INTERPOLATED_INVALID.
1676: If plex->interpolated == DMPLEX_INTERPOLATED_INVALID, DMPlexIsInterpolated_Internal() is called.
1677: It checks the actual topology and sets plex->interpolated on each rank separately to one of
1678: DMPLEX_INTERPOLATED_NONE, DMPLEX_INTERPOLATED_PARTIAL or DMPLEX_INTERPOLATED_FULL.
1680: If plex->interpolated != DMPLEX_INTERPOLATED_INVALID, this function just returns plex->interpolated.
1682: DMPlexInterpolate() sets plex->interpolated = DMPLEX_INTERPOLATED_FULL,
1683: and DMPlexUninterpolate() sets plex->interpolated = DMPLEX_INTERPOLATED_NONE.
1685: .seealso: DMPlexInterpolate(), DMPlexIsInterpolatedCollective()
1686: @*/
1687: PetscErrorCode DMPlexIsInterpolated(DM dm, DMPlexInterpolatedFlag *interpolated)
1688: {
1689: DM_Plex *plex = (DM_Plex *) dm->data;
1690: PetscErrorCode ierr;
1695: if (plex->interpolated < 0) {
1696: DMPlexIsInterpolated_Internal(dm, &plex->interpolated);
1697: } else if (PetscDefined (USE_DEBUG)) {
1698: DMPlexInterpolatedFlag flg;
1700: DMPlexIsInterpolated_Internal(dm, &flg);
1701: if (flg != plex->interpolated) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Stashed DMPlexInterpolatedFlag %s is inconsistent with current %s", DMPlexInterpolatedFlags[plex->interpolated], DMPlexInterpolatedFlags[flg]);
1702: }
1703: *interpolated = plex->interpolated;
1704: return(0);
1705: }
1707: /*@
1708: DMPlexIsInterpolatedCollective - Find out to what extent the DMPlex is topologically interpolated (in collective manner).
1710: Collective
1712: Input Parameter:
1713: . dm - The DM object
1715: Output Parameter:
1716: . interpolated - Flag whether the DM is interpolated
1718: Level: intermediate
1720: Notes:
1721: Unlike DMPlexIsInterpolated(), this is collective so the results are guaranteed to be the same on all ranks.
1723: This function will return DMPLEX_INTERPOLATED_MIXED if the results of DMPlexIsInterpolated() are different on different ranks.
1725: Developer Notes:
1726: Initially, plex->interpolatedCollective = DMPLEX_INTERPOLATED_INVALID.
1728: If plex->interpolatedCollective == DMPLEX_INTERPOLATED_INVALID, this function calls DMPlexIsInterpolated() which sets plex->interpolated.
1729: MPI_Allreduce() is then called and collectively consistent flag plex->interpolatedCollective is set and returned;
1730: if plex->interpolated varies on different ranks, plex->interpolatedCollective = DMPLEX_INTERPOLATED_MIXED,
1731: otherwise sets plex->interpolatedCollective = plex->interpolated.
1733: If plex->interpolatedCollective != DMPLEX_INTERPOLATED_INVALID, this function just returns plex->interpolatedCollective.
1735: .seealso: DMPlexInterpolate(), DMPlexIsInterpolated()
1736: @*/
1737: PetscErrorCode DMPlexIsInterpolatedCollective(DM dm, DMPlexInterpolatedFlag *interpolated)
1738: {
1739: DM_Plex *plex = (DM_Plex *) dm->data;
1740: PetscBool debug=PETSC_FALSE;
1741: PetscErrorCode ierr;
1746: PetscOptionsGetBool(((PetscObject) dm)->options, ((PetscObject) dm)->prefix, "-dm_plex_is_interpolated_collective_debug", &debug, NULL);
1747: if (plex->interpolatedCollective < 0) {
1748: DMPlexInterpolatedFlag min, max;
1749: MPI_Comm comm;
1751: PetscObjectGetComm((PetscObject)dm, &comm);
1752: DMPlexIsInterpolated(dm, &plex->interpolatedCollective);
1753: MPI_Allreduce(&plex->interpolatedCollective, &min, 1, MPIU_ENUM, MPI_MIN, comm);
1754: MPI_Allreduce(&plex->interpolatedCollective, &max, 1, MPIU_ENUM, MPI_MAX, comm);
1755: if (min != max) plex->interpolatedCollective = DMPLEX_INTERPOLATED_MIXED;
1756: if (debug) {
1757: PetscMPIInt rank;
1759: MPI_Comm_rank(comm, &rank);
1760: PetscSynchronizedPrintf(comm, "[%d] interpolated=%s interpolatedCollective=%s\n", rank, DMPlexInterpolatedFlags[plex->interpolated], DMPlexInterpolatedFlags[plex->interpolatedCollective]);
1761: PetscSynchronizedFlush(comm, PETSC_STDOUT);
1762: }
1763: }
1764: *interpolated = plex->interpolatedCollective;
1765: return(0);
1766: }