Actual source code: plexrefregular.c
1: #include <petsc/private/dmplextransformimpl.h>
3: /* Regular Refinment of Hybrid Meshes
5: We would like to express regular refinement as a small set of rules that can be applied on every point of the Plex
6: to automatically generate a refined Plex. In fact, we would like these rules to be general enough to encompass other
7: transformations, such as changing from one type of cell to another, as simplex to hex.
9: To start, we can create a function that takes an original cell type and returns the number of new cells replacing it
10: and the types of the new cells.
12: We need the group multiplication table for group actions from the dihedral group for each cell type.
14: We need an operator which takes in a cell, and produces a new set of cells with new faces and correct orientations. I think
15: we can just write this operator for faces with identity, and then compose the face orientation actions to get the actual
16: (face, orient) pairs for each subcell.
17: */
19: /*@
20: DMPlexRefineRegularGetAffineFaceTransforms - Gets the affine map from the reference face cell to each face in the given cell
22: Input Parameters:
23: + cr - The DMPlexCellRefiner object
24: - ct - The cell type
26: Output Parameters:
27: + Nf - The number of faces for this cell type
28: . v0 - The translation of the first vertex for each face
29: . J - The Jacobian for each face (map from original cell to subcell)
30: . invJ - The inverse Jacobian for each face
31: - detJ - The determinant of the Jacobian for each face
33: Note: The Jacobian and inverse Jacboian will be rectangular, and the inverse is really a generalized inverse.
34: $ v0 + j x_face = x_cell
35: $ invj (x_cell - v0) = x_face
37: Level: developer
39: .seealso: DMPlexCellRefinerGetAffineTransforms(), Create()
40: @*/
41: PetscErrorCode DMPlexRefineRegularGetAffineFaceTransforms(DMPlexTransform tr, DMPolytopeType ct, PetscInt *Nf, PetscReal *v0[], PetscReal *J[], PetscReal *invJ[], PetscReal *detJ[])
42: {
43: /*
44: 2
45: |\
46: | \
47: | \
48: | \
49: | \
50: | \
51: | \
52: 2 1
53: | \
54: | \
55: | \
56: 0---0-------1
57: v0[Nf][dc]: 3 x 2
58: J[Nf][df][dc]: 3 x 1 x 2
59: invJ[Nf][dc][df]: 3 x 2 x 1
60: detJ[Nf]: 3
61: */
62: static PetscReal tri_v0[] = {0.0, -1.0, 0.0, 0.0, -1.0, 0.0};
63: static PetscReal tri_J[] = {1.0, 0.0, -1.0, 1.0, 0.0, -1.0};
64: static PetscReal tri_invJ[] = {1.0, 0.0, -0.5, 0.5, 0.0, -1.0};
65: static PetscReal tri_detJ[] = {1.0, 1.414213562373095, 1.0};
66: /*
67: 3---------2---------2
68: | |
69: | |
70: | |
71: 3 1
72: | |
73: | |
74: | |
75: 0---------0---------1
77: v0[Nf][dc]: 4 x 2
78: J[Nf][df][dc]: 4 x 1 x 2
79: invJ[Nf][dc][df]: 4 x 2 x 1
80: detJ[Nf]: 4
81: */
82: static PetscReal quad_v0[] = {0.0, -1.0, 1.0, 0.0, 0.0, 1.0 -1.0, 0.0};
83: static PetscReal quad_J[] = {1.0, 0.0, 0.0, 1.0, -1.0, 0.0, 0.0, -1.0};
84: static PetscReal quad_invJ[] = {1.0, 0.0, 0.0, 1.0, -1.0, 0.0, 0.0, -1.0};
85: static PetscReal quad_detJ[] = {1.0, 1.0, 1.0, 1.0};
87: switch (ct) {
88: case DM_POLYTOPE_TRIANGLE: if (Nf) *Nf = 3; if (v0) *v0 = tri_v0; if (J) *J = tri_J; if (invJ) *invJ = tri_invJ; if (detJ) *detJ = tri_detJ; break;
89: case DM_POLYTOPE_QUADRILATERAL: if (Nf) *Nf = 4; if (v0) *v0 = quad_v0; if (J) *J = quad_J; if (invJ) *invJ = quad_invJ; if (detJ) *detJ = quad_detJ; break;
90: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported polytope type %s", DMPolytopeTypes[ct]);
91: }
92: return 0;
93: }
95: /*@
96: DMPlexRefineRegularGetAffineTransforms - Gets the affine map from the reference cell to each subcell
98: Input Parameters:
99: + cr - The DMPlexCellRefiner object
100: - ct - The cell type
102: Output Parameters:
103: + Nc - The number of subcells produced from this cell type
104: . v0 - The translation of the first vertex for each subcell
105: . J - The Jacobian for each subcell (map from reference cell to subcell)
106: - invJ - The inverse Jacobian for each subcell
108: Level: developer
110: .seealso: DMPlexRefineRegularGetAffineFaceTransforms(), DMPLEXREFINEREGULAR
111: @*/
112: PetscErrorCode DMPlexRefineRegularGetAffineTransforms(DMPlexTransform tr, DMPolytopeType ct, PetscInt *Nc, PetscReal *v0[], PetscReal *J[], PetscReal *invJ[])
113: {
114: /*
115: 2
116: |\
117: | \
118: | \
119: | \
120: | C \
121: | \
122: | \
123: 2---1---1
124: |\ D / \
125: | 2 0 \
126: |A \ / B \
127: 0---0-------1
128: */
129: static PetscReal tri_v0[] = {-1.0, -1.0, 0.0, -1.0, -1.0, 0.0, 0.0, -1.0};
130: static PetscReal tri_J[] = {0.5, 0.0,
131: 0.0, 0.5,
133: 0.5, 0.0,
134: 0.0, 0.5,
136: 0.5, 0.0,
137: 0.0, 0.5,
139: 0.0, -0.5,
140: 0.5, 0.5};
141: static PetscReal tri_invJ[] = {2.0, 0.0,
142: 0.0, 2.0,
144: 2.0, 0.0,
145: 0.0, 2.0,
147: 2.0, 0.0,
148: 0.0, 2.0,
150: 2.0, 2.0,
151: -2.0, 0.0};
152: /*
153: 3---------2---------2
154: | | |
155: | D 2 C |
156: | | |
157: 3----3----0----1----1
158: | | |
159: | A 0 B |
160: | | |
161: 0---------0---------1
162: */
163: static PetscReal quad_v0[] = {-1.0, -1.0, 0.0, -1.0, 0.0, 0.0, -1.0, 0.0};
164: static PetscReal quad_J[] = {0.5, 0.0,
165: 0.0, 0.5,
167: 0.5, 0.0,
168: 0.0, 0.5,
170: 0.5, 0.0,
171: 0.0, 0.5,
173: 0.5, 0.0,
174: 0.0, 0.5};
175: static PetscReal quad_invJ[] = {2.0, 0.0,
176: 0.0, 2.0,
178: 2.0, 0.0,
179: 0.0, 2.0,
181: 2.0, 0.0,
182: 0.0, 2.0,
184: 2.0, 0.0,
185: 0.0, 2.0};
186: /*
187: Bottom (viewed from top) Top
188: 1---------2---------2 7---------2---------6
189: | | | | | |
190: | B 2 C | | H 2 G |
191: | | | | | |
192: 3----3----0----1----1 3----3----0----1----1
193: | | | | | |
194: | A 0 D | | E 0 F |
195: | | | | | |
196: 0---------0---------3 4---------0---------5
197: */
198: static PetscReal hex_v0[] = {-1.0, -1.0, -1.0, -1.0, 0.0, -1.0, 0.0, 0.0, -1.0, 0.0, -1.0, -1.0,
199: -1.0, -1.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, -1.0, 0.0, 0.0};
200: static PetscReal hex_J[] = {0.5, 0.0, 0.0,
201: 0.0, 0.5, 0.0,
202: 0.0, 0.0, 0.5,
204: 0.5, 0.0, 0.0,
205: 0.0, 0.5, 0.0,
206: 0.0, 0.0, 0.5,
208: 0.5, 0.0, 0.0,
209: 0.0, 0.5, 0.0,
210: 0.0, 0.0, 0.5,
212: 0.5, 0.0, 0.0,
213: 0.0, 0.5, 0.0,
214: 0.0, 0.0, 0.5,
216: 0.5, 0.0, 0.0,
217: 0.0, 0.5, 0.0,
218: 0.0, 0.0, 0.5,
220: 0.5, 0.0, 0.0,
221: 0.0, 0.5, 0.0,
222: 0.0, 0.0, 0.5,
224: 0.5, 0.0, 0.0,
225: 0.0, 0.5, 0.0,
226: 0.0, 0.0, 0.5,
228: 0.5, 0.0, 0.0,
229: 0.0, 0.5, 0.0,
230: 0.0, 0.0, 0.5};
231: static PetscReal hex_invJ[] = {2.0, 0.0, 0.0,
232: 0.0, 2.0, 0.0,
233: 0.0, 0.0, 2.0,
235: 2.0, 0.0, 0.0,
236: 0.0, 2.0, 0.0,
237: 0.0, 0.0, 2.0,
239: 2.0, 0.0, 0.0,
240: 0.0, 2.0, 0.0,
241: 0.0, 0.0, 2.0,
243: 2.0, 0.0, 0.0,
244: 0.0, 2.0, 0.0,
245: 0.0, 0.0, 2.0,
247: 2.0, 0.0, 0.0,
248: 0.0, 2.0, 0.0,
249: 0.0, 0.0, 2.0,
251: 2.0, 0.0, 0.0,
252: 0.0, 2.0, 0.0,
253: 0.0, 0.0, 2.0,
255: 2.0, 0.0, 0.0,
256: 0.0, 2.0, 0.0,
257: 0.0, 0.0, 2.0,
259: 2.0, 0.0, 0.0,
260: 0.0, 2.0, 0.0,
261: 0.0, 0.0, 2.0};
263: switch (ct) {
264: case DM_POLYTOPE_TRIANGLE: if (Nc) *Nc = 4; if (v0) *v0 = tri_v0; if (J) *J = tri_J; if (invJ) *invJ = tri_invJ; break;
265: case DM_POLYTOPE_QUADRILATERAL: if (Nc) *Nc = 4; if (v0) *v0 = quad_v0; if (J) *J = quad_J; if (invJ) *invJ = quad_invJ; break;
266: case DM_POLYTOPE_HEXAHEDRON: if (Nc) *Nc = 8; if (v0) *v0 = hex_v0; if (J) *J = hex_J; if (invJ) *invJ = hex_invJ; break;
267: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported polytope type %s", DMPolytopeTypes[ct]);
268: }
269: return 0;
270: }
272: #if 0
273: static PetscErrorCode DMPlexCellRefinerGetCellVertices_Regular(DMPlexCellRefiner cr, DMPolytopeType ct, PetscInt *Nv, PetscReal *subcellV[])
274: {
275: static PetscReal seg_v[] = {-1.0, 0.0, 1.0};
276: static PetscReal tri_v[] = {-1.0, -1.0, 1.0, -1.0, -1.0, 1.0, 0.0, -1.0, 0.0, 0.0, -1.0, 0.0};
277: static PetscReal quad_v[] = {-1.0, -1.0, 1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 0.0, -1.0, 1.0, 0.0, 0.0, 1.0, -1.0, 0.0, 0.0, 0.0};
278: static PetscReal tet_v[] = {-1.0, -1.0, -1.0, 0.0, -1.0, -1.0, 1.0, -1.0, -1.0,
279: -1.0, 0.0, -1.0, 0.0, 0.0, -1.0, -1.0, 1.0, -1.0,
280: -1.0, -1.0, 0.0, 0.0, -1.0, 0.0, -1.0, 0.0, 0.0, -1.0, -1.0, 1.0};
281: static PetscReal hex_v[] = {-1.0, -1.0, -1.0, 0.0, -1.0, -1.0, 1.0, -1.0, -1.0,
282: -1.0, 0.0, -1.0, 0.0, 0.0, -1.0, 1.0, 0.0, -1.0,
283: -1.0, 1.0, -1.0, 0.0, 1.0, -1.0, 1.0, 1.0, -1.0,
284: -1.0, -1.0, 0.0, 0.0, -1.0, 0.0, 1.0, -1.0, 0.0,
285: -1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.0, 0.0,
286: -1.0, 1.0, 0.0, 0.0, 1.0, 0.0, 1.0, 1.0, 0.0,
287: -1.0, -1.0, 1.0, 0.0, -1.0, 1.0, 1.0, -1.0, 1.0,
288: -1.0, 0.0, 1.0, 0.0, 0.0, 1.0, 1.0, 0.0, 1.0,
289: -1.0, 1.0, 1.0, 0.0, 1.0, 1.0, 1.0, 1.0, 1.0};
291: switch (ct) {
292: case DM_POLYTOPE_SEGMENT: *Nv = 3; *subcellV = seg_v; break;
293: case DM_POLYTOPE_TRIANGLE: *Nv = 6; *subcellV = tri_v; break;
294: case DM_POLYTOPE_QUADRILATERAL: *Nv = 9; *subcellV = quad_v; break;
295: case DM_POLYTOPE_TETRAHEDRON: *Nv = 10; *subcellV = tet_v; break;
296: case DM_POLYTOPE_HEXAHEDRON: *Nv = 27; *subcellV = hex_v; break;
297: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No subcell vertices for cell type %s", DMPolytopeTypes[ct]);
298: }
299: return 0;
300: }
302: static PetscErrorCode DMPlexCellRefinerGetSubcellVertices_Regular(DMPlexCellRefiner cr, DMPolytopeType ct, DMPolytopeType rct, PetscInt r, PetscInt *Nv, PetscInt *subcellV[])
303: {
304: static PetscInt seg_v[] = {0, 1, 1, 2};
305: static PetscInt tri_v[] = {0, 3, 5, 3, 1, 4, 5, 4, 2, 3, 4, 5};
306: static PetscInt quad_v[] = {0, 4, 8, 7, 4, 1, 5, 8, 8, 5, 2, 6, 7, 8, 6, 3};
307: static PetscInt tet_v[] = {0, 3, 1, 6, 3, 2, 4, 8, 1, 4, 5, 7, 6, 8, 7, 9,
308: 1, 6, 3, 7, 8, 4, 3, 7, 7, 3, 1, 4, 7, 3, 8, 6};
309: static PetscInt hex_v[] = {0, 3, 4, 1, 9, 10, 13, 12, 3, 6, 7, 4, 12, 13, 16, 15, 4, 7, 8, 5, 13, 14, 17, 16, 1, 4 , 5 , 2, 10, 11, 14, 13,
310: 9, 12, 13, 10, 18, 19, 22, 21, 10, 13, 14, 11, 19, 20, 23, 22, 13, 16, 17, 14, 22, 23, 26, 25, 12, 15, 16, 13, 21, 22, 25, 24};
313: switch (ct) {
314: case DM_POLYTOPE_SEGMENT: *Nv = 2; *subcellV = &seg_v[r*(*Nv)]; break;
315: case DM_POLYTOPE_TRIANGLE: *Nv = 3; *subcellV = &tri_v[r*(*Nv)]; break;
316: case DM_POLYTOPE_QUADRILATERAL: *Nv = 4; *subcellV = &quad_v[r*(*Nv)]; break;
317: case DM_POLYTOPE_TETRAHEDRON: *Nv = 4; *subcellV = &tet_v[r*(*Nv)]; break;
318: case DM_POLYTOPE_HEXAHEDRON: *Nv = 8; *subcellV = &hex_v[r*(*Nv)]; break;
319: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No subcell vertices for cell type %s", DMPolytopeTypes[ct]);
320: }
321: return 0;
322: }
323: #endif
325: static PetscErrorCode DMPlexTransformView_Regular(DMPlexTransform tr, PetscViewer viewer)
326: {
327: PetscBool isascii;
331: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
332: if (isascii) {
333: const char *name;
335: PetscObjectGetName((PetscObject) tr, &name);
336: PetscViewerASCIIPrintf(viewer, "Regular refinement %s\n", name ? name : "");
337: } else {
338: SETERRQ(PetscObjectComm((PetscObject) tr), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMPlexTransform writing", ((PetscObject) viewer)->type_name);
339: }
340: return 0;
341: }
343: static PetscErrorCode DMPlexTransformSetUp_Regular(DMPlexTransform tr)
344: {
345: return 0;
346: }
348: static PetscErrorCode DMPlexTransformDestroy_Regular(DMPlexTransform tr)
349: {
350: DMPlexRefine_Regular *f = (DMPlexRefine_Regular *) tr->data;
352: PetscFree(f);
353: return 0;
354: }
356: PetscErrorCode DMPlexTransformGetSubcellOrientation_Regular(DMPlexTransform tr, DMPolytopeType sct, PetscInt sp, PetscInt so, DMPolytopeType tct, PetscInt r, PetscInt o, PetscInt *rnew, PetscInt *onew)
357: {
358: static PetscInt seg_seg[] = {1, -1, 0, -1,
359: 0, 0, 1, 0};
360: static PetscInt tri_seg[] = {2, -1, 1, -1, 0, -1,
361: 1, -1, 0, -1, 2, -1,
362: 0, -1, 2, -1, 1, -1,
363: 0, 0, 1, 0, 2, 0,
364: 1, 0, 2, 0, 0, 0,
365: 2, 0, 0, 0, 1, 0};
366: static PetscInt tri_tri[] = {1, -3, 0, -3, 2, -3, 3, -2,
367: 0, -2, 2, -2, 1, -2, 3, -1,
368: 2, -1, 1, -1, 0, -1, 3, -3,
369: 0, 0, 1, 0, 2, 0, 3, 0,
370: 1, 1, 2, 1, 0, 1, 3, 1,
371: 2, 2, 0, 2, 1, 2, 3, 2};
372: static PetscInt quad_seg[] = {1, 0, 0, 0, 3, 0, 2, 0,
373: 0, 0, 3, 0, 2, 0, 1, 0,
374: 3, 0, 2, 0, 1, 0, 0, 0,
375: 2, 0, 1, 0, 0, 0, 3, 0,
376: 0, 0, 1, 0, 2, 0, 3, 0,
377: 1, 0, 2, 0, 3, 0, 0, 0,
378: 2, 0, 3, 0, 0, 0, 1, 0,
379: 3, 0, 0, 0, 1, 0, 2, 0};
380: static PetscInt quad_quad[] = {2, -4, 1, -4, 0, -4, 3, -4,
381: 1, -3, 0, -3, 3, -3, 2, -3,
382: 0, -2, 3, -2, 2, -2, 1, -2,
383: 3, -1, 2, -1, 1, -1, 0, -1,
384: 0, 0, 1, 0, 2, 0, 3, 0,
385: 1, 1, 2, 1, 3, 1, 0, 1,
386: 2, 2, 3, 2, 0, 2, 1, 2,
387: 3, 3, 0, 3, 1, 3, 2, 3};
388: static PetscInt tseg_seg[] = {0, -1,
389: 0, 0,
390: 0, 0,
391: 0, -1};
392: static PetscInt tseg_tseg[] = {1, -2, 0, -2,
393: 1, -1, 0, -1,
394: 0, 0, 1, 0,
395: 0, 1, 1, 1};
396: static PetscInt tet_seg[] = {0, -1,
397: 0, 0,
398: 0, 0,
399: 0, -1,
400: 0, 0,
401: 0, 0,
402: 0, 0,
403: 0, 0,
404: 0, 0,
405: 0, 0,
406: 0, 0,
407: 0, 0,
408: 0, 0,
409: 0, 0,
410: 0, 0,
411: 0, 0,
412: 0, 0,
413: 0, 0,
414: 0, -1,
415: 0, 0,
416: 0, 0,
417: 0, -1,
418: 0, 0,
419: 0, 0};
420: static PetscInt tet_tri[] = {2, -1, 3, -1, 1, -3, 0, -2, 6, 1, 7, -3, 5, 2, 4, -3,
421: 3, -1, 1, -1, 2, -3, 0, -1, 4, 0, 5, 0, 6, 0, 7, 0,
422: 1, -1, 2, -1, 3, -3, 0, -3, 4, 0, 5, 0, 6, 0, 7, 0,
423: 3, -2, 2, -3, 0, -1, 1, -1, 7, -3, 6, 1, 4, 2, 5, -3,
424: 2, -3, 0, -2, 3, -1, 1, -3, 4, 0, 5, 0, 6, 0, 7, 0,
425: 0, -2, 3, -2, 2, -2, 1, -2, 4, 0, 5, 0, 6, 0, 7, 0,
426: 0, -1, 1, -2, 3, -2, 2, -2, 7, 1, 6, -3, 5, -3, 4, 2,
427: 1, -2, 3, -3, 0, -3, 2, -1, 4, 0, 5, 0, 6, 0, 7, 0,
428: 3, -3, 0, -1, 1, -1, 2, -3, 4, 0, 5, 0, 6, 0, 7, 0,
429: 1, -3, 0, -3, 2, -1, 3, -3, 6, -3, 7, 1, 4, -3, 5, 2,
430: 0, -3, 2, -2, 1, -2, 3, -2, 4, 0, 5, 0, 6, 0, 7, 0,
431: 2, -2, 1, -3, 0, -2, 3, -1, 4, 0, 5, 0, 6, 0, 7, 0,
432: 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, 6, 0, 7, 0,
433: 1, 0, 2, 2, 0, 1, 3, 1, 4, 0, 5, 0, 6, 0, 7, 0,
434: 2, 2, 0, 0, 1, 1, 3, 2, 4, 0, 5, 0, 6, 0, 7, 0,
435: 1, 2, 0, 1, 3, 1, 2, 2, 5, 0, 4, 0, 7, -1, 6, -1,
436: 0, 1, 3, 0, 1, 0, 2, 0, 4, 0, 5, 0, 6, 0, 7, 0,
437: 3, 0, 1, 2, 0, 2, 2, 1, 4, 0, 5, 0, 6, 0, 7, 0,
438: 2, 0, 3, 2, 0, 0, 1, 1, 4, -2, 5, -2, 7, 0, 6, 0,
439: 3, 2, 0, 2, 2, 1, 1, 2, 4, 0, 5, 0, 6, 0, 7, 0,
440: 0, 2, 2, 0, 3, 0, 1, 0, 4, 0, 5, 0, 6, 0, 7, 0,
441: 3, 1, 2, 1, 1, 2, 0, 2, 5, -2, 4, -2, 6, -1, 7, -1,
442: 2, 1, 1, 1, 3, 2, 0, 0, 4, 0, 5, 0, 6, 0, 7, 0,
443: 1, 1, 3, 1, 2, 2, 0, 1, 4, 0, 5, 0, 6, 0, 7, 0};
444: static PetscInt tet_tet[] = {2, -12, 3, -12, 1, -12, 0, -12, 6, -9, 7, -9, 5, -12, 4, -12,
445: 3, -11, 1, -11, 2, -11, 0, -11, 4, 0, 5, 0, 6, 0, 7, 0,
446: 1, -10, 2, -10, 3, -10, 0, -10, 4, 0, 5, 0, 6, 0, 7, 0,
447: 3, -9, 2, -9, 0, -9, 1, -9, 7, -9, 6, -9, 4, -12, 5, -12,
448: 2, -8, 0, -8, 3, -8, 1, -8, 4, 0, 5, 0, 6, 0, 7, 0,
449: 0, -7, 3, -7, 2, -7, 1, -7, 4, 0, 5, 0, 6, 0, 7, 0,
450: 0, -6, 1, -6, 3, -6, 2, -6, 4, -3, 5, -3, 7, -6, 6, -6,
451: 1, -5, 3, -5, 0, -5, 2, -5, 4, 0, 5, 0, 6, 0, 7, 0,
452: 3, -4, 0, -4, 1, -4, 2, -4, 4, 0, 5, 0, 6, 0, 7, 0,
453: 1, -3, 0, -3, 2, -3, 3, -3, 5, -3, 4, -3, 6, -6, 7, -6,
454: 0, -2, 2, -2, 1, -2, 3, -2, 4, 0, 5, 0, 6, 0, 7, 0,
455: 2, -1, 1, -1, 0, -1, 3, -1, 4, 0, 5, 0, 6, 0, 7, 0,
456: 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, 6, 0, 7, 0,
457: 1, 1, 2, 1, 0, 1, 3, 1, 4, 0, 5, 0, 6, 0, 7, 0,
458: 2, 2, 0, 2, 1, 2, 3, 2, 4, 0, 5, 0, 6, 0, 7, 0,
459: 1, 3, 0, 3, 3, 3, 2, 3, 5, 0, 4, 0, 7, 0, 6, 0,
460: 0, 4, 3, 4, 1, 4, 2, 4, 4, 0, 5, 0, 6, 0, 7, 0,
461: 3, 5, 1, 5, 0, 5, 2, 5, 4, 0, 5, 0, 6, 0, 7, 0,
462: 2, 6, 3, 6, 0, 6, 1, 6, 6, 6, 7, 6, 4, 6, 5, 6,
463: 3, 7, 0, 7, 2, 7, 1, 7, 4, 0, 5, 0, 6, 0, 7, 0,
464: 0, 8, 2, 8, 3, 8, 1, 8, 4, 0, 5, 0, 6, 0, 7, 0,
465: 3, 9, 2, 9, 1, 9, 0, 9, 7, 6, 6, 6, 5, 6, 4, 6,
466: 2, 10, 1, 10, 3, 10, 0, 10, 4, 0, 5, 0, 6, 0, 7, 0,
467: 1, 11, 3, 11, 2, 11, 0, 11, 4, 0, 5, 0, 6, 0, 7, 0};
468: static PetscInt hex_seg[] = {2, 0, 3, 0, 4, 0, 5, 0, 1, 0, 0, 0,
469: 4, 0, 5, 0, 0, 0, 1, 0, 3, 0, 2, 0,
470: 5, 0, 4, 0, 1, 0, 0, 0, 3, 0, 2, 0,
471: 3, 0, 2, 0, 4, 0, 5, 0, 0, 0, 1, 0,
472: 3, 0, 2, 0, 5, 0, 4, 0, 1, 0, 0, 0,
473: 4, 0, 5, 0, 1, 0, 0, 0, 2, 0, 3, 0,
474: 2, 0, 3, 0, 5, 0, 4, 0, 0, 0, 1, 0,
475: 5, 0, 4, 0, 0, 0, 1, 0, 2, 0, 3, 0,
476: 4, 0, 5, 0, 3, 0, 2, 0, 1, 0, 0, 0,
477: 5, 0, 4, 0, 3, 0, 2, 0, 0, 0, 1, 0,
478: 3, 0, 2, 0, 1, 0, 0, 0, 4, 0, 5, 0,
479: 2, 0, 3, 0, 0, 0, 1, 0, 4, 0, 5, 0,
480: 1, 0, 0, 0, 4, 0, 5, 0, 3, 0, 2, 0,
481: 1, 0, 0, 0, 5, 0, 4, 0, 2, 0, 3, 0,
482: 5, 0, 4, 0, 2, 0, 3, 0, 1, 0, 0, 0,
483: 1, 0, 0, 0, 2, 0, 3, 0, 4, 0, 5, 0,
484: 4, 0, 5, 0, 2, 0, 3, 0, 0, 0, 1, 0,
485: 3, 0, 2, 0, 0, 0, 1, 0, 5, 0, 4, 0,
486: 1, 0, 0, 0, 3, 0, 2, 0, 5, 0, 4, 0,
487: 2, 0, 3, 0, 1, 0, 0, 0, 5, 0, 4, 0,
488: 0, 0, 1, 0, 4, 0, 5, 0, 2, 0, 3, 0,
489: 0, 0, 1, 0, 3, 0, 2, 0, 4, 0, 5, 0,
490: 0, 0, 1, 0, 5, 0, 4, 0, 3, 0, 2, 0,
491: 0, 0, 1, 0, 2, 0, 3, 0, 5, 0, 4, 0,
492: 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0,
493: 0, 0, 1, 0, 5, 0, 4, 0, 2, 0, 3, 0,
494: 0, 0, 1, 0, 3, 0, 2, 0, 5, 0, 4, 0,
495: 0, 0, 1, 0, 4, 0, 5, 0, 3, 0, 2, 0,
496: 2, 0, 3, 0, 1, 0, 0, 0, 4, 0, 5, 0,
497: 1, 0, 0, 0, 3, 0, 2, 0, 4, 0, 5, 0,
498: 3, 0, 2, 0, 0, 0, 1, 0, 4, 0, 5, 0,
499: 4, 0, 5, 0, 2, 0, 3, 0, 1, 0, 0, 0,
500: 1, 0, 0, 0, 2, 0, 3, 0, 5, 0, 4, 0,
501: 5, 0, 4, 0, 2, 0, 3, 0, 0, 0, 1, 0,
502: 1, 0, 0, 0, 5, 0, 4, 0, 3, 0, 2, 0,
503: 1, 0, 0, 0, 4, 0, 5, 0, 2, 0, 3, 0,
504: 2, 0, 3, 0, 0, 0, 1, 0, 5, 0, 4, 0,
505: 3, 0, 2, 0, 1, 0, 0, 0, 5, 0, 4, 0,
506: 5, 0, 4, 0, 3, 0, 2, 0, 1, 0, 0, 0,
507: 4, 0, 5, 0, 3, 0, 2, 0, 0, 0, 1, 0,
508: 5, 0, 4, 0, 0, 0, 1, 0, 3, 0, 2, 0,
509: 2, 0, 3, 0, 5, 0, 4, 0, 1, 0, 0, 0,
510: 4, 0, 5, 0, 1, 0, 0, 0, 3, 0, 2, 0,
511: 3, 0, 2, 0, 5, 0, 4, 0, 0, 0, 1, 0,
512: 3, 0, 2, 0, 4, 0, 5, 0, 1, 0, 0, 0,
513: 5, 0, 4, 0, 1, 0, 0, 0, 2, 0, 3, 0,
514: 4, 0, 5, 0, 0, 0, 1, 0, 2, 0, 3, 0,
515: 2, 0, 3, 0, 4, 0, 5, 0, 0, 0, 1, 0};
516: static PetscInt hex_quad[] = {7, -2, 4, -2, 5, -2, 6, -2, 8, -3, 11, -3, 10, -3, 9, -3, 3, 1, 2, 1, 1, 1, 0, 1,
517: 8, -2, 9, -2, 10, -2, 11, -2, 3, -4, 0, -4, 1, -4, 2, -4, 7, 0, 4, 0, 5, 0, 6, 0,
518: 9, 1, 8, 1, 11, 1, 10, 1, 0, 3, 3, 3, 2, 3, 1, 3, 5, 2, 6, 2, 7, 2, 4, 2,
519: 6, 3, 5, 3, 4, 3, 7, 3, 10, -1, 9, -1, 8, -1, 11, -1, 2, -4, 3, -4, 0, -4, 1, -4,
520: 4, 1, 7, 1, 6, 1, 5, 1, 11, 2, 8, 2, 9, 2, 10, 2, 1, 3, 0, 3, 3, 3, 2, 3,
521: 10, -4, 11, -4, 8, -4, 9, -4, 2, 1, 1, 1, 0, 1, 3, 1, 6, -1, 5, -1, 4, -1, 7, -1,
522: 5, -4, 6, -4, 7, -4, 4, -4, 9, 0, 10, 0, 11, 0, 8, 0, 0, -2, 1, -2, 2, -2, 3, -2,
523: 11, 3, 10, 3, 9, 3, 8, 3, 1, -2, 2, -2, 3, -2, 0, -2, 4, -3, 7, -3, 6, -3, 5, -3,
524: 11, -1, 8, -1, 9, -1, 10, -1, 7, 3, 4, 3, 5, 3, 6, 3, 2, 2, 1, 2, 0, 2, 3, 2,
525: 10, 2, 9, 2, 8, 2, 11, 2, 5, 1, 6, 1, 7, 1, 4, 1, 1, -1, 2, -1, 3, -1, 0, -1,
526: 5, 2, 4, 2, 7, 2, 6, 2, 1, 2, 0, 2, 3, 2, 2, 2, 10, -4, 9, -4, 8, -4, 11, -4,
527: 4, -3, 5, -3, 6, -3, 7, -3, 0, -3, 1, -3, 2, -3, 3, -3, 8, -2, 11, -2, 10, -2, 9, -2,
528: 3, 1, 0, 1, 1, 1, 2, 1, 9, -4, 8, -4, 11, -4, 10, -4, 6, 3, 7, 3, 4, 3, 5, 3,
529: 1, 3, 2, 3, 3, 3, 0, 3, 10, 1, 11, 1, 8, 1, 9, 1, 5, -4, 4, -4, 7, -4, 6, -4,
530: 8, 0, 11, 0, 10, 0, 9, 0, 4, -4, 7, -4, 6, -4, 5, -4, 0, 0, 3, 0, 2, 0, 1, 0,
531: 0, 0, 1, 0, 2, 0, 3, 0, 5, -1, 4, -1, 7, -1, 6, -1, 9, -3, 8, -3, 11, -3, 10, -3,
532: 9, -3, 10, -3, 11, -3, 8, -3, 6, -2, 5, -2, 4, -2, 7, -2, 3, -3, 0, -3, 1, -3, 2, -3,
533: 7, 0, 6, 0, 5, 0, 4, 0, 2, -1, 3, -1, 0, -1, 1, -1, 11, 3, 8, 3, 9, 3, 10, 3,
534: 2, 2, 3, 2, 0, 2, 1, 2, 6, 2, 7, 2, 4, 2, 5, 2, 10, 2, 11, 2, 8, 2, 9, 2,
535: 6, -1, 7, -1, 4, -1, 5, -1, 3, 0, 2, 0, 1, 0, 0, 0, 9, 1, 10, 1, 11, 1, 8, 1,
536: 2, -4, 1, -4, 0, -4, 3, -4, 11, -2, 10, -2, 9, -2, 8, -2, 7, -2, 6, -2, 5, -2, 4, -2,
537: 1, -1, 0, -1, 3, -1, 2, -1, 4, 0, 5, 0, 6, 0, 7, 0, 11, -1, 10, -1, 9, -1, 8, -1,
538: 0, -2, 3, -2, 2, -2, 1, -2, 8, 3, 9, 3, 10, 3, 11, 3, 4, 1, 5, 1, 6, 1, 7, 1,
539: 3, -3, 2, -3, 1, -3, 0, -3, 7, -3, 6, -3, 5, -3, 4, -3, 8, 0, 9, 0, 10, 0, 11, 0,
540: 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, 6, 0, 7, 0, 8, 0, 9, 0, 10, 0, 11, 0,
541: 1, 3, 2, 3, 3, 3, 0, 3, 11, -2, 10, -2, 9, -2, 8, -2, 4, 1, 5, 1, 6, 1, 7, 1,
542: 2, 2, 3, 2, 0, 2, 1, 2, 7, -3, 6, -3, 5, -3, 4, -3, 11, -1, 10, -1, 9, -1, 8, -1,
543: 3, 1, 0, 1, 1, 1, 2, 1, 8, 3, 9, 3, 10, 3, 11, 3, 7, -2, 6, -2, 5, -2, 4, -2,
544: 5, 2, 4, 2, 7, 2, 6, 2, 0, -3, 1, -3, 2, -3, 3, -3, 9, 1, 10, 1, 11, 1, 8, 1,
545: 1, -1, 0, -1, 3, -1, 2, -1, 5, -1, 4, -1, 7, -1, 6, -1, 10, 2, 11, 2, 8, 2, 9, 2,
546: 4, -3, 5, -3, 6, -3, 7, -3, 1, 2, 0, 2, 3, 2, 2, 2, 11, 3, 8, 3, 9, 3, 10, 3,
547: 8, 0, 11, 0, 10, 0, 9, 0, 7, 3, 4, 3, 5, 3, 6, 3, 3, -3, 0, -3, 1, -3, 2, -3,
548: 3, -3, 2, -3, 1, -3, 0, -3, 6, 2, 7, 2, 4, 2, 5, 2, 9, -3, 8, -3, 11, -3, 10, -3,
549: 9, -3, 10, -3, 11, -3, 8, -3, 5, 1, 6, 1, 7, 1, 4, 1, 0, 0, 3, 0, 2, 0, 1, 0,
550: 0, -2, 3, -2, 2, -2, 1, -2, 9, -4, 8, -4, 11, -4, 10, -4, 5, -4, 4, -4, 7, -4, 6, -4,
551: 2, -4, 1, -4, 0, -4, 3, -4, 10, 1, 11, 1, 8, 1, 9, 1, 6, 3, 7, 3, 4, 3, 5, 3,
552: 7, 0, 6, 0, 5, 0, 4, 0, 3, 0, 2, 0, 1, 0, 0, 0, 8, -2, 11, -2, 10, -2, 9, -2,
553: 6, -1, 7, -1, 4, -1, 5, -1, 2, -1, 3, -1, 0, -1, 1, -1, 10, -4, 9, -4, 8, -4, 11, -4,
554: 11, -1, 8, -1, 9, -1, 10, -1, 4, -4, 7, -4, 6, -4, 5, -4, 1, -1, 2, -1, 3, -1, 0, -1,
555: 10, 2, 9, 2, 8, 2, 11, 2, 6, -2, 5, -2, 4, -2, 7, -2, 2, 2, 1, 2, 0, 2, 3, 2,
556: 8, -2, 9, -2, 10, -2, 11, -2, 0, 3, 3, 3, 2, 3, 1, 3, 4, -3, 7, -3, 6, -3, 5, -3,
557: 4, 1, 7, 1, 6, 1, 5, 1, 8, -3, 11, -3, 10, -3, 9, -3, 0, -2, 1, -2, 2, -2, 3, -2,
558: 9, 1, 8, 1, 11, 1, 10, 1, 3, -4, 0, -4, 1, -4, 2, -4, 6, -1, 5, -1, 4, -1, 7, -1,
559: 5, -4, 6, -4, 7, -4, 4, -4, 10, -1, 9, -1, 8, -1, 11, -1, 1, 3, 0, 3, 3, 3, 2, 3,
560: 7, -2, 4, -2, 5, -2, 6, -2, 11, 2, 8, 2, 9, 2, 10, 2, 2, -4, 3, -4, 0, -4, 1, -4,
561: 10, -4, 11, -4, 8, -4, 9, -4, 1, -2, 2, -2, 3, -2, 0, -2, 5, 2, 6, 2, 7, 2, 4, 2,
562: 11, 3, 10, 3, 9, 3, 8, 3, 2, 1, 1, 1, 0, 1, 3, 1, 7, 0, 4, 0, 5, 0, 6, 0,
563: 6, 3, 5, 3, 4, 3, 7, 3, 9, 0, 10, 0, 11, 0, 8, 0, 3, 1, 2, 1, 1, 1, 0, 1};
564: static PetscInt hex_hex[] = {3, -24, 0, -24, 4, -24, 5, -24, 2, -24, 6, -24, 7, -24, 1, -24,
565: 3, -23, 5, -23, 6, -23, 2, -23, 0, -23, 1, -23, 7, -23, 4, -23,
566: 4, -22, 0, -22, 1, -22, 7, -22, 5, -22, 6, -22, 2, -22, 3, -22,
567: 6, -21, 7, -21, 1, -21, 2, -21, 5, -21, 3, -21, 0, -21, 4, -21,
568: 1, -20, 2, -20, 6, -20, 7, -20, 0, -20, 4, -20, 5, -20, 3, -20,
569: 6, -19, 2, -19, 3, -19, 5, -19, 7, -19, 4, -19, 0, -19, 1, -19,
570: 4, -18, 5, -18, 3, -18, 0, -18, 7, -18, 1, -18, 2, -18, 6, -18,
571: 1, -17, 7, -17, 4, -17, 0, -17, 2, -17, 3, -17, 5, -17, 6, -17,
572: 2, -16, 3, -16, 5, -16, 6, -16, 1, -16, 7, -16, 4, -16, 0, -16,
573: 7, -15, 4, -15, 0, -15, 1, -15, 6, -15, 2, -15, 3, -15, 5, -15,
574: 7, -14, 1, -14, 2, -14, 6, -14, 4, -14, 5, -14, 3, -14, 0, -14,
575: 0, -13, 4, -13, 5, -13, 3, -13, 1, -13, 2, -13, 6, -13, 7, -13,
576: 5, -12, 4, -12, 7, -12, 6, -12, 3, -12, 2, -12, 1, -12, 0, -12,
577: 7, -11, 6, -11, 5, -11, 4, -11, 1, -11, 0, -11, 3, -11, 2, -11,
578: 0, -10, 1, -10, 7, -10, 4, -10, 3, -10, 5, -10, 6, -10, 2, -10,
579: 4, -9, 7, -9, 6, -9, 5, -9, 0, -9, 3, -9, 2, -9, 1, -9,
580: 5, -8, 6, -8, 2, -8, 3, -8, 4, -8, 0, -8, 1, -8, 7, -8,
581: 2, -7, 6, -7, 7, -7, 1, -7, 3, -7, 0, -7, 4, -7, 5, -7,
582: 6, -6, 5, -6, 4, -6, 7, -6, 2, -6, 1, -6, 0, -6, 3, -6,
583: 5, -5, 3, -5, 0, -5, 4, -5, 6, -5, 7, -5, 1, -5, 2, -5,
584: 2, -4, 1, -4, 0, -4, 3, -4, 6, -4, 5, -4, 4, -4, 7, -4,
585: 1, -3, 0, -3, 3, -3, 2, -3, 7, -3, 6, -3, 5, -3, 4, -3,
586: 0, -2, 3, -2, 2, -2, 1, -2, 4, -2, 7, -2, 6, -2, 5, -2,
587: 3, -1, 2, -1, 1, -1, 0, -1, 5, -1, 4, -1, 7, -1, 6, -1,
588: 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, 6, 0, 7, 0,
589: 1, 1, 2, 1, 3, 1, 0, 1, 7, 1, 4, 1, 5, 1, 6, 1,
590: 2, 2, 3, 2, 0, 2, 1, 2, 6, 2, 7, 2, 4, 2, 5, 2,
591: 3, 3, 0, 3, 1, 3, 2, 3, 5, 3, 6, 3, 7, 3, 4, 3,
592: 4, 4, 0, 4, 3, 4, 5, 4, 7, 4, 6, 4, 2, 4, 1, 4,
593: 7, 5, 4, 5, 5, 5, 6, 5, 1, 5, 2, 5, 3, 5, 0, 5,
594: 1, 6, 7, 6, 6, 6, 2, 6, 0, 6, 3, 6, 5, 6, 4, 6,
595: 3, 7, 2, 7, 6, 7, 5, 7, 0, 7, 4, 7, 7, 7, 1, 7,
596: 5, 8, 6, 8, 7, 8, 4, 8, 3, 8, 0, 8, 1, 8, 2, 8,
597: 4, 9, 7, 9, 1, 9, 0, 9, 5, 9, 3, 9, 2, 9, 6, 9,
598: 4, 10, 5, 10, 6, 10, 7, 10, 0, 10, 1, 10, 2, 10, 3, 10,
599: 6, 11, 7, 11, 4, 11, 5, 11, 2, 11, 3, 11, 0, 11, 1, 11,
600: 3, 12, 5, 12, 4, 12, 0, 12, 2, 12, 1, 12, 7, 12, 6, 12,
601: 6, 13, 2, 13, 1, 13, 7, 13, 5, 13, 4, 13, 0, 13, 3, 13,
602: 1, 14, 0, 14, 4, 14, 7, 14, 2, 14, 6, 14, 5, 14, 3, 14,
603: 6, 15, 5, 15, 3, 15, 2, 15, 7, 15, 1, 15, 0, 15, 4, 15,
604: 0, 16, 4, 16, 7, 16, 1, 16, 3, 16, 2, 16, 6, 16, 5, 16,
605: 0, 17, 3, 17, 5, 17, 4, 17, 1, 17, 7, 17, 6, 17, 2, 17,
606: 5, 18, 3, 18, 2, 18, 6, 18, 4, 18, 7, 18, 1, 18, 0, 18,
607: 7, 19, 6, 19, 2, 19, 1, 19, 4, 19, 0, 19, 3, 19, 5, 19,
608: 2, 20, 1, 20, 7, 20, 6, 20, 3, 20, 5, 20, 4, 20, 0, 20,
609: 7, 21, 1, 21, 0, 21, 4, 21, 6, 21, 5, 21, 3, 21, 2, 21,
610: 2, 22, 6, 22, 5, 22, 3, 22, 1, 22, 0, 22, 4, 22, 7, 22,
611: 5, 23, 4, 23, 0, 23, 3, 23, 6, 23, 2, 23, 1, 23, 7, 23};
612: static PetscInt trip_seg[] = {1, 0, 2, 0, 0, 0,
613: 2, 0, 0, 0, 1, 0,
614: 0, 0, 1, 0, 2, 0,
615: 0, -1, 2, -1, 1, -1,
616: 1, -1, 0, -1, 2, -1,
617: 2, -1, 1, -1, 0, -1,
618: 0, 0, 1, 0, 2, 0,
619: 2, 0, 0, 0, 1, 0,
620: 1, 0, 2, 0, 0, 0,
621: 2, -1, 1, -1, 0, -1,
622: 1, -1, 0, -1, 2, -1,
623: 0, -1, 2, -1, 1, -1};
624: static PetscInt trip_tri[] = {1, 1, 2, 1, 0, 1, 3, 1,
625: 2, 2, 0, 2, 1, 2, 3, 2,
626: 0, 0, 1, 0, 2, 0, 3, 0,
627: 2, -1, 1, -1, 0, -1, 3, -3,
628: 0, -2, 2, -2, 1, -2, 3, -1,
629: 1, -3, 0, -3, 2, -3, 3, -2,
630: 0, 0, 1, 0, 2, 0, 3, 0,
631: 2, 2, 0, 2, 1, 2, 3, 2,
632: 1, 1, 2, 1, 0, 1, 3, 1,
633: 1, -3, 0, -3, 2, -3, 3, -2,
634: 0, -2, 2, -2, 1, -2, 3, -1,
635: 2, -1, 1, -1, 0, -1, 3, -3};
636: static PetscInt trip_quad[] = {4, -1, 5, -1, 3, -1, 1, -1, 2, -1, 0, -1,
637: 5, -1, 3, -1, 4, -1, 2, -1, 0, -1, 1, -1,
638: 3, -1, 4, -1, 5, -1, 0, -1, 1, -1, 2, -1,
639: 0, -3, 2, -3, 1, -3, 3, -3, 5, -3, 4, -3,
640: 1, -3, 0, -3, 2, -3, 4, -3, 3, -3, 5, -3,
641: 2, -3, 1, -3, 0, -3, 5, -3, 4, -3, 3, -3,
642: 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0,
643: 2, 0, 0, 0, 1, 0, 5, 0, 3, 0, 4, 0,
644: 1, 0, 2, 0, 0, 0, 4, 0, 5, 0, 3, 0,
645: 5, 2, 4, 2, 3, 2, 2, 2, 1, 2, 0, 2,
646: 4, 2, 3, 2, 5, 2, 1, 2, 0, 2, 2, 2,
647: 3, 2, 5, 2, 4, 2, 0, 2, 2, 2, 1, 2};
648: static PetscInt trip_trip[] = {5, -6, 6, -6, 4, -6, 7, -6, 1, -6, 2, -6, 0, -6, 3, -6,
649: 6, -5, 4, -5, 5, -5, 7, -5, 2, -5, 0, -5, 1, -5, 3, -5,
650: 4, -4, 5, -4, 6, -4, 7, -4, 0, -4, 1, -4, 2, -4, 3, -4,
651: 2, -3, 1, -3, 0, -3, 3, -1, 6, -3, 5, -3, 4, -3, 7, -1,
652: 0, -2, 2, -2, 1, -2, 3, -3, 4, -2, 6, -2, 5, -2, 7, -3,
653: 1, -1, 0, -1, 2, -1, 3, -2, 5, -1, 4, -1, 6, -1, 7, -2,
654: 0, 0, 1, 0, 2, 0, 3, 0, 4, 0, 5, 0, 6, 0, 7, 0,
655: 2, 1, 0, 1, 1, 1, 3, 1, 6, 1, 4, 1, 5, 1, 7, 1,
656: 1, 2, 2, 2, 0, 2, 3, 2, 5, 2, 6, 2, 4, 2, 7, 2,
657: 5, 3, 4, 3, 6, 3, 7, 4, 1, 3, 0, 3, 2, 3, 3, 4,
658: 4, 4, 6, 4, 5, 4, 7, 5, 0, 4, 2, 4, 1, 4, 3, 5,
659: 6, 5, 5, 5, 4, 5, 7, 3, 2, 5, 1, 5, 0, 5, 3, 3};
660: static PetscInt ttri_tseg[] = {2, -2, 1, -2, 0, -2,
661: 1, -2, 0, -2, 2, -2,
662: 0, -2, 2, -2, 1, -2,
663: 2, -1, 1, -1, 0, -1,
664: 1, -1, 0, -1, 2, -1,
665: 0, -1, 2, -1, 1, -1,
666: 0, 0, 1, 0, 2, 0,
667: 1, 0, 2, 0, 0, 0,
668: 2, 0, 0, 0, 1, 0,
669: 0, 1, 1, 1, 2, 1,
670: 1, 1, 2, 1, 0, 1,
671: 2, 1, 0, 1, 1, 1};
672: static PetscInt ttri_ttri[] = {1, -6, 0, -6, 2, -6, 3, -5,
673: 0, -5, 2, -5, 1, -5, 3, -4,
674: 2, -4, 1, -4, 0, -4, 3, -6,
675: 1, -3, 0, -3, 2, -3, 3, -2,
676: 0, -2, 2, -2, 1, -2, 3, -1,
677: 2, -1, 1, -1, 0, -1, 3, -3,
678: 0, 0, 1, 0, 2, 0, 3, 0,
679: 1, 1, 2, 1, 0, 1, 3, 1,
680: 2, 2, 0, 2, 1, 2, 3, 2,
681: 0, 3, 1, 3, 2, 3, 3, 3,
682: 1, 4, 2, 4, 0, 4, 3, 4,
683: 2, 5, 0, 5, 1, 5, 3, 5};
684: static PetscInt tquad_tvert[] = {0, -1,
685: 0, -1,
686: 0, -1,
687: 0, -1,
688: 0, 0,
689: 0, 0,
690: 0, 0,
691: 0, 0,
692: 0, 0,
693: 0, 0,
694: 0, 0,
695: 0, 0,
696: 0, -1,
697: 0, -1,
698: 0, -1,
699: 0, -1};
700: static PetscInt tquad_tseg[] = {1, 1, 0, 1, 3, 1, 2, 1,
701: 0, 1, 3, 1, 2, 1, 1, 1,
702: 3, 1, 2, 1, 1, 1, 0, 1,
703: 2, 1, 1, 1, 0, 1, 3, 1,
704: 1, 0, 0, 0, 3, 0, 2, 0,
705: 0, 0, 3, 0, 2, 0, 1, 0,
706: 3, 0, 2, 0, 1, 0, 0, 0,
707: 2, 0, 1, 0, 0, 0, 3, 0,
708: 0, 0, 1, 0, 2, 0, 3, 0,
709: 1, 0, 2, 0, 3, 0, 0, 0,
710: 2, 0, 3, 0, 0, 0, 1, 0,
711: 3, 0, 0, 0, 1, 0, 2, 0,
712: 0, 1, 1, 1, 2, 1, 3, 1,
713: 1, 1, 2, 1, 3, 1, 0, 1,
714: 2, 1, 3, 1, 0, 1, 1, 1,
715: 3, 1, 0, 1, 1, 1, 2, 1};
716: static PetscInt tquad_tquad[] = {2, -8, 1, -8, 0, -8, 3, -8,
717: 1, -7, 0, -7, 3, -7, 2, -7,
718: 0, -6, 3, -6, 2, -6, 1, -6,
719: 3, -5, 2, -5, 1, -5, 0, -5,
720: 2, -4, 1, -4, 0, -4, 3, -4,
721: 1, -3, 0, -3, 3, -3, 2, -3,
722: 0, -2, 3, -2, 2, -2, 1, -2,
723: 3, -1, 2, -1, 1, -1, 0, -1,
724: 0, 0, 1, 0, 2, 0, 3, 0,
725: 1, 1, 2, 1, 3, 1, 0, 1,
726: 2, 2, 3, 2, 0, 2, 1, 2,
727: 3, 3, 0, 3, 1, 3, 2, 3,
728: 0, 4, 1, 4, 2, 4, 3, 4,
729: 1, 5, 2, 5, 3, 5, 0, 5,
730: 2, 6, 3, 6, 0, 6, 1, 6,
731: 3, 7, 0, 7, 1, 7, 2, 7};
734: *rnew = r; *onew = o;
735: if (!so) return 0;
736: switch (sct) {
737: case DM_POLYTOPE_POINT: break;
738: case DM_POLYTOPE_POINT_PRISM_TENSOR:
739: *onew = so < 0 ? -(o+1) : o;
740: break;
741: case DM_POLYTOPE_SEGMENT:
742: switch (tct) {
743: case DM_POLYTOPE_POINT: break;
744: case DM_POLYTOPE_SEGMENT:
745: *rnew = seg_seg[(so+1)*4 + r*2];
746: *onew = DMPolytopeTypeComposeOrientation(tct, o, seg_seg[(so+1)*4 + r*2 + 1]);
747: break;
748: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
749: }
750: break;
751: case DM_POLYTOPE_TRIANGLE:
752: switch (tct) {
753: case DM_POLYTOPE_SEGMENT:
754: *rnew = tri_seg[(so+3)*6 + r*2];
755: *onew = DMPolytopeTypeComposeOrientation(tct, o, tri_seg[(so+3)*6 + r*2 + 1]);
756: break;
757: case DM_POLYTOPE_TRIANGLE:
758: *rnew = tri_tri[(so+3)*8 + r*2];
759: *onew = DMPolytopeTypeComposeOrientation(tct, o, tri_tri[(so+3)*8 + r*2 + 1]);
760: break;
761: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
762: }
763: break;
764: case DM_POLYTOPE_QUADRILATERAL:
765: switch (tct) {
766: case DM_POLYTOPE_POINT: break;
767: case DM_POLYTOPE_SEGMENT:
768: *rnew = quad_seg[(so+4)*8 + r*2];
769: *onew = DMPolytopeTypeComposeOrientation(tct, o, quad_seg[(so+4)*8 + r*2 + 1]);
770: break;
771: case DM_POLYTOPE_QUADRILATERAL:
772: *rnew = quad_quad[(so+4)*8 + r*2];
773: *onew = DMPolytopeTypeComposeOrientation(tct, o, quad_quad[(so+4)*8 + r*2 + 1]);
774: break;
775: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
776: }
777: break;
778: case DM_POLYTOPE_SEG_PRISM_TENSOR:
779: switch (tct) {
780: case DM_POLYTOPE_POINT_PRISM_TENSOR:
781: *rnew = tseg_seg[(so+2)*2 + r*2];
782: *onew = DMPolytopeTypeComposeOrientation(tct, o, tseg_seg[(so+2)*2 + r*2 + 1]);
783: break;
784: case DM_POLYTOPE_SEG_PRISM_TENSOR:
785: *rnew = tseg_tseg[(so+2)*4 + r*2];
786: *onew = DMPolytopeTypeComposeOrientation(tct, o, tseg_tseg[(so+2)*4 + r*2 + 1]);
787: break;
788: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
789: }
790: break;
791: case DM_POLYTOPE_TETRAHEDRON:
792: switch (tct) {
793: case DM_POLYTOPE_SEGMENT:
794: *rnew = tet_seg[(so+12)*2 + r*2];
795: *onew = DMPolytopeTypeComposeOrientation(tct, o, tet_seg[(so+12)*2 + r*2 + 1]);
796: break;
797: case DM_POLYTOPE_TRIANGLE:
798: *rnew = tet_tri[(so+12)*16 + r*2];
799: *onew = DMPolytopeTypeComposeOrientation(tct, o, tet_tri[(so+12)*16 + r*2 + 1]);
800: break;
801: case DM_POLYTOPE_TETRAHEDRON:
802: *rnew = tet_tet[(so+12)*16 + r*2];
803: *onew = DMPolytopeTypeComposeOrientation(tct, o, tet_tet[(so+12)*16 + r*2 + 1]);
804: break;
805: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
806: }
807: break;
808: case DM_POLYTOPE_HEXAHEDRON:
809: switch (tct) {
810: case DM_POLYTOPE_POINT: break;
811: case DM_POLYTOPE_SEGMENT:
812: *rnew = hex_seg[(so+24)*12 + r*2];
813: *onew = DMPolytopeTypeComposeOrientation(tct, o, hex_seg[(so+24)*12 + r*2 + 1]);
814: break;
815: case DM_POLYTOPE_QUADRILATERAL:
816: *rnew = hex_quad[(so+24)*24 + r*2];
817: *onew = DMPolytopeTypeComposeOrientation(tct, o, hex_quad[(so+24)*24 + r*2 + 1]);
818: break;
819: case DM_POLYTOPE_HEXAHEDRON:
820: *rnew = hex_hex[(so+24)*16 + r*2];
821: *onew = DMPolytopeTypeComposeOrientation(tct, o, hex_hex[(so+24)*16 + r*2 + 1]);
822: break;
823: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
824: }
825: break;
826: case DM_POLYTOPE_TRI_PRISM:
827: switch (tct) {
828: case DM_POLYTOPE_SEGMENT:
829: *rnew = trip_seg[(so+6)*6 + r*2];
830: *onew = DMPolytopeTypeComposeOrientation(tct, o, trip_seg[(so+6)*6 + r*2 + 1]);
831: break;
832: case DM_POLYTOPE_TRIANGLE:
833: *rnew = trip_tri[(so+6)*8 + r*2];
834: *onew = DMPolytopeTypeComposeOrientation(tct, o, trip_tri[(so+6)*8 + r*2 + 1]);
835: break;
836: case DM_POLYTOPE_QUADRILATERAL:
837: *rnew = trip_quad[(so+6)*12 + r*2];
838: *onew = DMPolytopeTypeComposeOrientation(tct, o, trip_quad[(so+6)*12 + r*2 + 1]);
839: break;
840: case DM_POLYTOPE_TRI_PRISM:
841: *rnew = trip_trip[(so+6)*16 + r*2];
842: *onew = DMPolytopeTypeComposeOrientation(tct, o, trip_trip[(so+6)*16 + r*2 + 1]);
843: break;
844: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
845: }
846: break;
847: case DM_POLYTOPE_TRI_PRISM_TENSOR:
848: switch (tct) {
849: case DM_POLYTOPE_SEG_PRISM_TENSOR:
850: *rnew = ttri_tseg[(so+6)*6 + r*2];
851: *onew = DMPolytopeTypeComposeOrientation(tct, o, ttri_tseg[(so+6)*6 + r*2 + 1]);
852: break;
853: case DM_POLYTOPE_TRI_PRISM_TENSOR:
854: *rnew = ttri_ttri[(so+6)*8 + r*2];
855: *onew = DMPolytopeTypeComposeOrientation(tct, o, ttri_ttri[(so+6)*8 + r*2 + 1]);
856: break;
857: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
858: }
859: break;
860: case DM_POLYTOPE_QUAD_PRISM_TENSOR:
861: switch (tct) {
862: case DM_POLYTOPE_POINT_PRISM_TENSOR:
863: *rnew = tquad_tvert[(so+8)*2 + r*2];
864: *onew = DMPolytopeTypeComposeOrientation(tct, o, tquad_tvert[(so+8)*2 + r*2 + 1]);
865: break;
866: case DM_POLYTOPE_SEG_PRISM_TENSOR:
867: *rnew = tquad_tseg[(so+8)*8 + r*2];
868: *onew = DMPolytopeTypeComposeOrientation(tct, o, tquad_tseg[(so+8)*8 + r*2 + 1]);
869: break;
870: case DM_POLYTOPE_QUAD_PRISM_TENSOR:
871: *rnew = tquad_tquad[(so+8)*8 + r*2];
872: *onew = DMPolytopeTypeComposeOrientation(tct, o, tquad_tquad[(so+8)*8 + r*2 + 1]);
873: break;
874: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type %s is not produced by %s", DMPolytopeTypes[tct], DMPolytopeTypes[sct]);
875: }
876: break;
877: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell type %s", DMPolytopeTypes[sct]);
878: }
879: return 0;
880: }
882: PetscErrorCode DMPlexTransformCellRefine_Regular(DMPlexTransform tr, DMPolytopeType source, PetscInt p, PetscInt *rt, PetscInt *Nt, DMPolytopeType *target[], PetscInt *size[], PetscInt *cone[], PetscInt *ornt[])
883: {
884: /* All vertices remain in the refined mesh */
885: static DMPolytopeType vertexT[] = {DM_POLYTOPE_POINT};
886: static PetscInt vertexS[] = {1};
887: static PetscInt vertexC[] = {0};
888: static PetscInt vertexO[] = {0};
889: /* Split all edges with a new vertex, making two new 2 edges
890: 0--0--0--1--1
891: */
892: static DMPolytopeType segT[] = {DM_POLYTOPE_POINT, DM_POLYTOPE_SEGMENT};
893: static PetscInt segS[] = {1, 2};
894: static PetscInt segC[] = {DM_POLYTOPE_POINT, 1, 0, 0, DM_POLYTOPE_POINT, 0, 0, DM_POLYTOPE_POINT, 0, 0, DM_POLYTOPE_POINT, 1, 1, 0};
895: static PetscInt segO[] = { 0, 0, 0, 0};
896: /* Do not split tensor edges */
897: static DMPolytopeType tvertT[] = {DM_POLYTOPE_POINT_PRISM_TENSOR};
898: static PetscInt tvertS[] = {1};
899: static PetscInt tvertC[] = {DM_POLYTOPE_POINT, 1, 0, 0, DM_POLYTOPE_POINT, 1, 1, 0};
900: static PetscInt tvertO[] = { 0, 0};
901: /* Add 3 edges inside every triangle, making 4 new triangles.
902: 2
903: |\
904: | \
905: | \
906: 0 1
907: | C \
908: | \
909: | \
910: 2---1---1
911: |\ D / \
912: 1 2 0 0
913: |A \ / B \
914: 0-0-0---1---1
915: */
916: static DMPolytopeType triT[] = {DM_POLYTOPE_SEGMENT, DM_POLYTOPE_TRIANGLE};
917: static PetscInt triS[] = {3, 4};
918: static PetscInt triC[] = {DM_POLYTOPE_POINT, 1, 0, 0, DM_POLYTOPE_POINT, 1, 1, 0,
919: DM_POLYTOPE_POINT, 1, 1, 0, DM_POLYTOPE_POINT, 1, 2, 0,
920: DM_POLYTOPE_POINT, 1, 2, 0, DM_POLYTOPE_POINT, 1, 0, 0,
921: DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 0, 2, DM_POLYTOPE_SEGMENT, 1, 2, 1,
922: DM_POLYTOPE_SEGMENT, 1, 0, 1, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_SEGMENT, 0, 0,
923: DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 1, 1, 1, DM_POLYTOPE_SEGMENT, 1, 2, 0,
924: DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 0, 2};
925: static PetscInt triO[] = {0, 0,
926: 0, 0,
927: 0, 0,
928: 0, -1, 0,
929: 0, 0, -1,
930: -1, 0, 0,
931: 0, 0, 0};
932: /* Add a vertex in the center of each quadrilateral, and 4 edges inside, making 4 new quads.
933: 3----1----2----0----2
934: | | |
935: 0 D 2 C 1
936: | | |
937: 3----3----0----1----1
938: | | |
939: 1 A 0 B 0
940: | | |
941: 0----0----0----1----1
942: */
943: static DMPolytopeType quadT[] = {DM_POLYTOPE_POINT, DM_POLYTOPE_SEGMENT, DM_POLYTOPE_QUADRILATERAL};
944: static PetscInt quadS[] = {1, 4, 4};
945: static PetscInt quadC[] = {DM_POLYTOPE_POINT, 1, 0, 0, DM_POLYTOPE_POINT, 0, 0,
946: DM_POLYTOPE_POINT, 1, 1, 0, DM_POLYTOPE_POINT, 0, 0,
947: DM_POLYTOPE_POINT, 1, 2, 0, DM_POLYTOPE_POINT, 0, 0,
948: DM_POLYTOPE_POINT, 1, 3, 0, DM_POLYTOPE_POINT, 0, 0,
949: DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 0, 3, DM_POLYTOPE_SEGMENT, 1, 3, 1,
950: DM_POLYTOPE_SEGMENT, 1, 0, 1, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 0, 0,
951: DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 1, 1, 1, DM_POLYTOPE_SEGMENT, 1, 2, 0, DM_POLYTOPE_SEGMENT, 0, 2,
952: DM_POLYTOPE_SEGMENT, 0, 3, DM_POLYTOPE_SEGMENT, 0, 2, DM_POLYTOPE_SEGMENT, 1, 2, 1, DM_POLYTOPE_SEGMENT, 1, 3, 0};
953: static PetscInt quadO[] = {0, 0,
954: 0, 0,
955: 0, 0,
956: 0, 0,
957: 0, 0, -1, 0,
958: 0, 0, 0, -1,
959: -1, 0, 0, 0,
960: 0, -1, 0, 0};
961: /* Add 1 edge inside every tensor quad, making 2 new tensor quads
962: 2----2----1----3----3
963: | | |
964: | | |
965: | | |
966: 4 A 6 B 5
967: | | |
968: | | |
969: | | |
970: 0----0----0----1----1
971: */
972: static DMPolytopeType tsegT[] = {DM_POLYTOPE_POINT_PRISM_TENSOR, DM_POLYTOPE_SEG_PRISM_TENSOR};
973: static PetscInt tsegS[] = {1, 2};
974: static PetscInt tsegC[] = {DM_POLYTOPE_POINT, 1, 0, 0, DM_POLYTOPE_POINT, 1, 1, 0,
975: DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 2, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 0, 0,
976: DM_POLYTOPE_SEGMENT, 1, 0, 1, DM_POLYTOPE_SEGMENT, 1, 1, 1, DM_POLYTOPE_POINT_PRISM_TENSOR, 0, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 3, 0};
977: static PetscInt tsegO[] = {0, 0,
978: 0, 0, 0, 0,
979: 0, 0, 0, 0};
980: /* Add 1 edge and 8 triangles inside every cell, making 8 new tets
981: The vertices of our reference tet are [(-1, -1, -1), (-1, 1, -1), (1, -1, -1), (-1, -1, 1)], which we call [v0, v1, v2, v3]. The first
982: three edges are [v0, v1], [v1, v2], [v2, v0] called e0, e1, and e2, and then three edges to the top point [v0, v3], [v1, v3], [v2, v3]
983: called e3, e4, and e5. The faces of a tet, given in DMPlexGetRawFaces_Internal() are
984: [v0, v1, v2], [v0, v3, v1], [v0, v2, v3], [v2, v1, v3]
985: The first four tets just cut off the corners, using the replica notation for new vertices,
986: [v0, (e0, 0), (e2, 0), (e3, 0)]
987: [(e0, 0), v1, (e1, 0), (e4, 0)]
988: [(e2, 0), (e1, 0), v2, (e5, 0)]
989: [(e3, 0), (e4, 0), (e5, 0), v3 ]
990: The next four tets match a vertex to the newly created faces from cutting off those first tets.
991: [(e2, 0), (e3, 0), (e0, 0), (e5, 0)]
992: [(e4, 0), (e1, 0), (e0, 0), (e5, 0)]
993: [(e5, 0), (e0, 0), (e2, 0), (e1, 0)]
994: [(e5, 0), (e0, 0), (e4, 0), (e3, 0)]
995: We can see that a new edge is introduced in the cell [(e0, 0), (e5, 0)] which we call (-1, 0). The first four faces created are
996: [(e2, 0), (e0, 0), (e3, 0)]
997: [(e0, 0), (e1, 0), (e4, 0)]
998: [(e2, 0), (e5, 0), (e1, 0)]
999: [(e3, 0), (e4, 0), (e5, 0)]
1000: The next four, from the second group of tets, are
1001: [(e2, 0), (e0, 0), (e5, 0)]
1002: [(e4, 0), (e0, 0), (e5, 0)]
1003: [(e0, 0), (e1, 0), (e5, 0)]
1004: [(e5, 0), (e3, 0), (e0, 0)]
1005: I could write a program to generate these orientations by comparing the faces from GetRawFaces() with my existing table.
1006: */
1007: static DMPolytopeType tetT[] = {DM_POLYTOPE_SEGMENT, DM_POLYTOPE_TRIANGLE, DM_POLYTOPE_TETRAHEDRON};
1008: static PetscInt tetS[] = {1, 8, 8};
1009: static PetscInt tetC[] = {DM_POLYTOPE_POINT, 2, 0, 0, 0, DM_POLYTOPE_POINT, 2, 2, 1, 0,
1010: DM_POLYTOPE_SEGMENT, 1, 0, 2, DM_POLYTOPE_SEGMENT, 1, 1, 2, DM_POLYTOPE_SEGMENT, 1, 2, 2,
1011: DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 3, 0, DM_POLYTOPE_SEGMENT, 1, 1, 1,
1012: DM_POLYTOPE_SEGMENT, 1, 2, 0, DM_POLYTOPE_SEGMENT, 1, 3, 2, DM_POLYTOPE_SEGMENT, 1, 0, 1,
1013: DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_SEGMENT, 1, 3, 1, DM_POLYTOPE_SEGMENT, 1, 2, 1,
1014: DM_POLYTOPE_SEGMENT, 1, 0, 2, DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 1, 2, 0,
1015: DM_POLYTOPE_SEGMENT, 1, 1, 1, DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 1, 3, 1,
1016: DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 3, 2, DM_POLYTOPE_SEGMENT, 0, 0,
1017: DM_POLYTOPE_SEGMENT, 1, 2, 1, DM_POLYTOPE_SEGMENT, 1, 1, 2, DM_POLYTOPE_SEGMENT, 0, 0,
1018: DM_POLYTOPE_TRIANGLE, 1, 0, 0, DM_POLYTOPE_TRIANGLE, 1, 1, 0, DM_POLYTOPE_TRIANGLE, 1, 2, 0, DM_POLYTOPE_TRIANGLE, 0, 0,
1019: DM_POLYTOPE_TRIANGLE, 1, 0, 1, DM_POLYTOPE_TRIANGLE, 1, 1, 2, DM_POLYTOPE_TRIANGLE, 0, 1, DM_POLYTOPE_TRIANGLE, 1, 3, 1,
1020: DM_POLYTOPE_TRIANGLE, 1, 0, 2, DM_POLYTOPE_TRIANGLE, 0, 2, DM_POLYTOPE_TRIANGLE, 1, 2, 1, DM_POLYTOPE_TRIANGLE, 1, 3, 0,
1021: DM_POLYTOPE_TRIANGLE, 0, 3, DM_POLYTOPE_TRIANGLE, 1, 1, 1, DM_POLYTOPE_TRIANGLE, 1, 2, 2, DM_POLYTOPE_TRIANGLE, 1, 3, 2,
1022: DM_POLYTOPE_TRIANGLE, 0, 0, DM_POLYTOPE_TRIANGLE, 1, 2, 3, DM_POLYTOPE_TRIANGLE, 0, 4, DM_POLYTOPE_TRIANGLE, 0, 7,
1023: DM_POLYTOPE_TRIANGLE, 0, 1, DM_POLYTOPE_TRIANGLE, 1, 3, 3, DM_POLYTOPE_TRIANGLE, 0, 5, DM_POLYTOPE_TRIANGLE, 0, 6,
1024: DM_POLYTOPE_TRIANGLE, 0, 4, DM_POLYTOPE_TRIANGLE, 0, 6, DM_POLYTOPE_TRIANGLE, 0, 2, DM_POLYTOPE_TRIANGLE, 1, 0, 3,
1025: DM_POLYTOPE_TRIANGLE, 0, 5, DM_POLYTOPE_TRIANGLE, 0, 7, DM_POLYTOPE_TRIANGLE, 0, 3, DM_POLYTOPE_TRIANGLE, 1, 1, 3};
1026: static PetscInt tetO[] = {0, 0,
1027: 0, 0, 0,
1028: 0, 0, 0,
1029: 0, 0, 0,
1030: 0, 0, 0,
1031: 0, 0, -1,
1032: 0, 0, -1,
1033: 0, -1, -1,
1034: 0, -1, 0,
1035: 0, 0, 0, 0,
1036: 0, 0, 0, 0,
1037: 0, 0, 0, 0,
1038: 0, 0, 0, 0,
1039: -2, 0, 0, -1,
1040: -1, 1, 0, 0,
1041: -1, -1, -3, 2,
1042: -1, 0, -1, 1};
1043: /* Add a vertex in the center of each cell, add 6 edges and 12 quads inside every cell, making 8 new hexes
1044: The vertices of our reference hex are (-1, -1, -1), (-1, 1, -1), (1, 1, -1), (1, -1, -1), (-1, -1, 1), (1, -1, 1), (1, 1, 1), (-1, 1, 1) which we call [v0, v1, v2, v3, v4, v5, v6, v7]. The fours edges around the bottom [v0, v1], [v1, v2], [v2, v3], [v3, v0] are [e0, e1, e2, e3], and likewise around the top [v4, v5], [v5, v6], [v6, v7], [v7, v4] are [e4, e5, e6, e7]. Finally [v0, v4], [v1, v7], [v2, v6], [v3, v5] are [e9, e10, e11, e8]. The faces of a hex, given in DMPlexGetRawFaces_Internal(), oriented with outward normals, are
1045: [v0, v1, v2, v3] f0 bottom
1046: [v4, v5, v6, v7] f1 top
1047: [v0, v3, v5, v4] f2 front
1048: [v2, v1, v7, v6] f3 back
1049: [v3, v2, v6, v5] f4 right
1050: [v0, v4, v7, v1] f5 left
1051: The eight hexes are divided into four on the bottom, and four on the top,
1052: [v0, (e0, 0), (f0, 0), (e3, 0), (e9, 0), (f2, 0), (c0, 0), (f5, 0)]
1053: [(e0, 0), v1, (e1, 0), (f0, 0), (f5, 0), (c0, 0), (f3, 0), (e10, 0)]
1054: [(f0, 0), (e1, 0), v2, (e2, 0), (c0, 0), (f4, 0), (e11, 0), (f3, 0)]
1055: [(e3, 0), (f0, 0), (e2, 0), v3, (f2, 0), (e8, 0), (f4, 0), (c0, 0)]
1056: [(e9, 0), (f5, 0), (c0, 0), (f2, 0), v4, (e4, 0), (f1, 0), (e7, 0)]
1057: [(f2, 0), (c0, 0), (f4, 0), (e8, 0), (e4, 0), v5, (e5, 0), (f1, 0)]
1058: [(c0, 0), (f3, 0), (e11, 0), (f4, 0), (f1, 0), (e5, 0), v6, (e6, 0)]
1059: [(f5, 0), (e10, 0), (f3, 0), (c0, 0), (e7, 0), (f1, 0), (e6, 0), v7]
1060: The 6 internal edges will go from the faces to the central vertex. The 12 internal faces can be divided into groups of 4 by the plane on which they sit. First the faces on the x-y plane are,
1061: [(e9, 0), (f2, 0), (c0, 0), (f5, 0)]
1062: [(f5, 0), (c0, 0), (f3, 0), (e10, 0)]
1063: [(c0, 0), (f4, 0), (e11, 0), (f3, 0)]
1064: [(f2, 0), (e8, 0), (f4, 0), (c0, 0)]
1065: and on the x-z plane,
1066: [(f0, 0), (e0, 0), (f5, 0), (c0, 0)]
1067: [(c0, 0), (f5, 0), (e7, 0), (f1, 0)]
1068: [(f4, 0), (c0, 0), (f1, 0), (e5, 0)]
1069: [(e2, 0), (f0, 0), (c0, 0), (f4, 0)]
1070: and on the y-z plane,
1071: [(e3, 0), (f2, 0), (c0, 0), (f0, 0)]
1072: [(f2, 0), (e4, 0), (f1, 0), (c0, 0)]
1073: [(c0, 0), (f1, 0), (e6, 0), (f3, 0)]
1074: [(f0, 0), (c0, 0), (f3, 0), (e1, 0)]
1075: */
1076: static DMPolytopeType hexT[] = {DM_POLYTOPE_POINT, DM_POLYTOPE_SEGMENT, DM_POLYTOPE_QUADRILATERAL, DM_POLYTOPE_HEXAHEDRON};
1077: static PetscInt hexS[] = {1, 6, 12, 8};
1078: static PetscInt hexC[] = {DM_POLYTOPE_POINT, 1, 0, 0, DM_POLYTOPE_POINT, 0, 0,
1079: DM_POLYTOPE_POINT, 1, 1, 0, DM_POLYTOPE_POINT, 0, 0,
1080: DM_POLYTOPE_POINT, 1, 2, 0, DM_POLYTOPE_POINT, 0, 0,
1081: DM_POLYTOPE_POINT, 1, 3, 0, DM_POLYTOPE_POINT, 0, 0,
1082: DM_POLYTOPE_POINT, 1, 4, 0, DM_POLYTOPE_POINT, 0, 0,
1083: DM_POLYTOPE_POINT, 1, 5, 0, DM_POLYTOPE_POINT, 0, 0,
1084: DM_POLYTOPE_SEGMENT, 1, 2, 3, DM_POLYTOPE_SEGMENT, 0, 2, DM_POLYTOPE_SEGMENT, 0, 5, DM_POLYTOPE_SEGMENT, 1, 5, 0,
1085: DM_POLYTOPE_SEGMENT, 0, 5, DM_POLYTOPE_SEGMENT, 0, 3, DM_POLYTOPE_SEGMENT, 1, 3, 1, DM_POLYTOPE_SEGMENT, 1, 5, 2,
1086: DM_POLYTOPE_SEGMENT, 0, 4, DM_POLYTOPE_SEGMENT, 1, 4, 1, DM_POLYTOPE_SEGMENT, 1, 3, 3, DM_POLYTOPE_SEGMENT, 0, 3,
1087: DM_POLYTOPE_SEGMENT, 1, 2, 1, DM_POLYTOPE_SEGMENT, 1, 4, 3, DM_POLYTOPE_SEGMENT, 0, 4, DM_POLYTOPE_SEGMENT, 0, 2,
1088: DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 5, 3, DM_POLYTOPE_SEGMENT, 0, 5, DM_POLYTOPE_SEGMENT, 0, 0,
1089: DM_POLYTOPE_SEGMENT, 0, 5, DM_POLYTOPE_SEGMENT, 1, 5, 1, DM_POLYTOPE_SEGMENT, 1, 1, 3, DM_POLYTOPE_SEGMENT, 0, 1,
1090: DM_POLYTOPE_SEGMENT, 0, 4, DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 1, 1, 1, DM_POLYTOPE_SEGMENT, 1, 4, 2,
1091: DM_POLYTOPE_SEGMENT, 1, 0, 2, DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 0, 4, DM_POLYTOPE_SEGMENT, 1, 4, 0,
1092: DM_POLYTOPE_SEGMENT, 1, 2, 0, DM_POLYTOPE_SEGMENT, 0, 2, DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 1, 0, 3,
1093: DM_POLYTOPE_SEGMENT, 1, 2, 2, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 0, 2,
1094: DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 1, 1, 2, DM_POLYTOPE_SEGMENT, 1, 3, 2, DM_POLYTOPE_SEGMENT, 0, 3,
1095: DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 0, 3, DM_POLYTOPE_SEGMENT, 1, 3, 0, DM_POLYTOPE_SEGMENT, 1, 0, 1,
1096: DM_POLYTOPE_QUADRILATERAL, 1, 0, 0, DM_POLYTOPE_QUADRILATERAL, 0, 0, DM_POLYTOPE_QUADRILATERAL, 1, 2, 0, DM_POLYTOPE_QUADRILATERAL, 0, 4, DM_POLYTOPE_QUADRILATERAL, 0, 8, DM_POLYTOPE_QUADRILATERAL, 1, 5, 0,
1097: DM_POLYTOPE_QUADRILATERAL, 1, 0, 1, DM_POLYTOPE_QUADRILATERAL, 0, 1, DM_POLYTOPE_QUADRILATERAL, 0, 4, DM_POLYTOPE_QUADRILATERAL, 1, 3, 1, DM_POLYTOPE_QUADRILATERAL, 0, 11, DM_POLYTOPE_QUADRILATERAL, 1, 5, 3,
1098: DM_POLYTOPE_QUADRILATERAL, 1, 0, 2, DM_POLYTOPE_QUADRILATERAL, 0, 2, DM_POLYTOPE_QUADRILATERAL, 0, 7, DM_POLYTOPE_QUADRILATERAL, 1, 3, 0, DM_POLYTOPE_QUADRILATERAL, 1, 4, 1, DM_POLYTOPE_QUADRILATERAL, 0, 11,
1099: DM_POLYTOPE_QUADRILATERAL, 1, 0, 3, DM_POLYTOPE_QUADRILATERAL, 0, 3, DM_POLYTOPE_QUADRILATERAL, 1, 2, 1, DM_POLYTOPE_QUADRILATERAL, 0, 7, DM_POLYTOPE_QUADRILATERAL, 1, 4, 0, DM_POLYTOPE_QUADRILATERAL, 0, 8,
1100: DM_POLYTOPE_QUADRILATERAL, 0, 0, DM_POLYTOPE_QUADRILATERAL, 1, 1, 0, DM_POLYTOPE_QUADRILATERAL, 1, 2, 3, DM_POLYTOPE_QUADRILATERAL, 0, 5, DM_POLYTOPE_QUADRILATERAL, 0, 9, DM_POLYTOPE_QUADRILATERAL, 1, 5, 1,
1101: DM_POLYTOPE_QUADRILATERAL, 0, 3, DM_POLYTOPE_QUADRILATERAL, 1, 1, 1, DM_POLYTOPE_QUADRILATERAL, 1, 2, 2, DM_POLYTOPE_QUADRILATERAL, 0, 6, DM_POLYTOPE_QUADRILATERAL, 1, 4, 3, DM_POLYTOPE_QUADRILATERAL, 0, 9,
1102: DM_POLYTOPE_QUADRILATERAL, 0, 2, DM_POLYTOPE_QUADRILATERAL, 1, 1, 2, DM_POLYTOPE_QUADRILATERAL, 0, 6, DM_POLYTOPE_QUADRILATERAL, 1, 3, 3, DM_POLYTOPE_QUADRILATERAL, 1, 4, 2, DM_POLYTOPE_QUADRILATERAL, 0, 10,
1103: DM_POLYTOPE_QUADRILATERAL, 0, 1, DM_POLYTOPE_QUADRILATERAL, 1, 1, 3, DM_POLYTOPE_QUADRILATERAL, 0, 5, DM_POLYTOPE_QUADRILATERAL, 1, 3, 2, DM_POLYTOPE_QUADRILATERAL, 0, 10, DM_POLYTOPE_QUADRILATERAL, 1, 5, 2};
1104: static PetscInt hexO[] = {0, 0,
1105: 0, 0,
1106: 0, 0,
1107: 0, 0,
1108: 0, 0,
1109: 0, 0,
1110: 0, 0, -1, -1,
1111: 0, -1, -1, 0,
1112: -1, -1, 0, 0,
1113: -1, 0, 0, -1,
1114: -1, 0, 0, -1,
1115: -1, -1, 0, 0,
1116: 0, -1, -1, 0,
1117: 0, 0, -1, -1,
1118: 0, 0, -1, -1,
1119: -1, 0, 0, -1,
1120: -1, -1, 0, 0,
1121: 0, -1, -1, 0,
1122: 0, 0, 0, 0, -2, 0,
1123: 0, 0, -3, 0, -2, 0,
1124: 0, 0, -3, 0, 0, 0,
1125: 0, 0, 0, 0, 0, 0,
1126: -2, 0, 0, 0, -2, 0,
1127: -2, 0, 0, 0, 0, 0,
1128: -2, 0, -3, 0, 0, 0,
1129: -2, 0, -3, 0, -2, 0};
1130: /* Add 3 quads inside every triangular prism, making 4 new prisms. */
1131: static DMPolytopeType tripT[] = {DM_POLYTOPE_SEGMENT, DM_POLYTOPE_TRIANGLE, DM_POLYTOPE_QUADRILATERAL, DM_POLYTOPE_TRI_PRISM};
1132: static PetscInt tripS[] = {3, 4, 6, 8};
1133: static PetscInt tripC[] = {DM_POLYTOPE_POINT, 1, 2, 0, DM_POLYTOPE_POINT, 1, 3, 0,
1134: DM_POLYTOPE_POINT, 1, 3, 0, DM_POLYTOPE_POINT, 1, 4, 0,
1135: DM_POLYTOPE_POINT, 1, 4, 0, DM_POLYTOPE_POINT, 1, 2, 0,
1136: DM_POLYTOPE_SEGMENT, 1, 2, 3, DM_POLYTOPE_SEGMENT, 0, 2, DM_POLYTOPE_SEGMENT, 1, 4, 1,
1137: DM_POLYTOPE_SEGMENT, 1, 2, 1, DM_POLYTOPE_SEGMENT, 1, 3, 3, DM_POLYTOPE_SEGMENT, 0, 0,
1138: DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 1, 3, 1, DM_POLYTOPE_SEGMENT, 1, 4, 3,
1139: DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 0, 2,
1140: DM_POLYTOPE_SEGMENT, 1, 0, 1, DM_POLYTOPE_SEGMENT, 1, 3, 0, DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 1, 2, 0,
1141: DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 4, 0, DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 1, 3, 0,
1142: DM_POLYTOPE_SEGMENT, 1, 0, 2, DM_POLYTOPE_SEGMENT, 1, 2, 0, DM_POLYTOPE_SEGMENT, 0, 2, DM_POLYTOPE_SEGMENT, 1, 4, 0,
1143: DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 1, 3, 2, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_SEGMENT, 1, 2, 2,
1144: DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 1, 4, 2, DM_POLYTOPE_SEGMENT, 1, 1, 1, DM_POLYTOPE_SEGMENT, 1, 3, 2,
1145: DM_POLYTOPE_SEGMENT, 0, 2, DM_POLYTOPE_SEGMENT, 1, 2, 2, DM_POLYTOPE_SEGMENT, 1, 1, 2, DM_POLYTOPE_SEGMENT, 1, 4, 2,
1146: DM_POLYTOPE_TRIANGLE, 1, 0, 0, DM_POLYTOPE_TRIANGLE, 0, 0, DM_POLYTOPE_QUADRILATERAL, 1, 2, 0, DM_POLYTOPE_QUADRILATERAL, 0, 2, DM_POLYTOPE_QUADRILATERAL, 1, 4, 1,
1147: DM_POLYTOPE_TRIANGLE, 1, 0, 2, DM_POLYTOPE_TRIANGLE, 0, 1, DM_POLYTOPE_QUADRILATERAL, 1, 2, 1, DM_POLYTOPE_QUADRILATERAL, 1, 3, 0, DM_POLYTOPE_QUADRILATERAL, 0, 0,
1148: DM_POLYTOPE_TRIANGLE, 1, 0, 1, DM_POLYTOPE_TRIANGLE, 0, 2, DM_POLYTOPE_QUADRILATERAL, 0, 1, DM_POLYTOPE_QUADRILATERAL, 1, 3, 1, DM_POLYTOPE_QUADRILATERAL, 1, 4, 0,
1149: DM_POLYTOPE_TRIANGLE, 1, 0, 3, DM_POLYTOPE_TRIANGLE, 0, 3, DM_POLYTOPE_QUADRILATERAL, 0, 0, DM_POLYTOPE_QUADRILATERAL, 0, 1, DM_POLYTOPE_QUADRILATERAL, 0, 2,
1150: DM_POLYTOPE_TRIANGLE, 0, 0, DM_POLYTOPE_TRIANGLE, 1, 1, 0, DM_POLYTOPE_QUADRILATERAL, 1, 2, 3, DM_POLYTOPE_QUADRILATERAL, 0, 5, DM_POLYTOPE_QUADRILATERAL, 1, 4, 2,
1151: DM_POLYTOPE_TRIANGLE, 0, 1, DM_POLYTOPE_TRIANGLE, 1, 1, 1, DM_POLYTOPE_QUADRILATERAL, 1, 2, 2, DM_POLYTOPE_QUADRILATERAL, 1, 3, 3, DM_POLYTOPE_QUADRILATERAL, 0, 3,
1152: DM_POLYTOPE_TRIANGLE, 0, 2, DM_POLYTOPE_TRIANGLE, 1, 1, 2, DM_POLYTOPE_QUADRILATERAL, 0, 4, DM_POLYTOPE_QUADRILATERAL, 1, 3, 2, DM_POLYTOPE_QUADRILATERAL, 1, 4, 3,
1153: DM_POLYTOPE_TRIANGLE, 0, 3, DM_POLYTOPE_TRIANGLE, 1, 1, 3, DM_POLYTOPE_QUADRILATERAL, 0, 3, DM_POLYTOPE_QUADRILATERAL, 0, 4, DM_POLYTOPE_QUADRILATERAL, 0, 5};
1154: static PetscInt tripO[] = {0, 0,
1155: 0, 0,
1156: 0, 0,
1157: 0, -1, -1,
1158: -1, 0, -1,
1159: -1, -1, 0,
1160: 0, 0, 0,
1161: -1, 0, -1, -1,
1162: -1, 0, -1, -1,
1163: -1, 0, -1, -1,
1164: 0, -1, -1, 0,
1165: 0, -1, -1, 0,
1166: 0, -1, -1, 0,
1167: 0, 0, 0, -3, 0,
1168: 0, 0, 0, 0, -3,
1169: 0, 0, -3, 0, 0,
1170: 2, 0, 0, 0, 0,
1171: -2, 0, 0, -3, 0,
1172: -2, 0, 0, 0, -3,
1173: -2, 0, -3, 0, 0,
1174: -2, 0, 0, 0, 0};
1175: /* Add 3 tensor quads inside every tensor triangular prism, making 4 new prisms.
1176: 2
1177: |\
1178: | \
1179: | \
1180: 0---1
1182: 2
1184: 0 1
1186: 2
1187: |\
1188: | \
1189: | \
1190: 0---1
1191: */
1192: static DMPolytopeType ttriT[] = {DM_POLYTOPE_SEG_PRISM_TENSOR, DM_POLYTOPE_TRI_PRISM_TENSOR};
1193: static PetscInt ttriS[] = {3, 4};
1194: static PetscInt ttriC[] = {DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 2, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 3, 0,
1195: DM_POLYTOPE_SEGMENT, 1, 0, 1, DM_POLYTOPE_SEGMENT, 1, 1, 1, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 3, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 4, 0,
1196: DM_POLYTOPE_SEGMENT, 1, 0, 2, DM_POLYTOPE_SEGMENT, 1, 1, 2, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 4, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 2, 0,
1197: DM_POLYTOPE_TRIANGLE, 1, 0, 0, DM_POLYTOPE_TRIANGLE, 1, 1, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 2, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 2, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 4, 1,
1198: DM_POLYTOPE_TRIANGLE, 1, 0, 1, DM_POLYTOPE_TRIANGLE, 1, 1, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 2, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 3, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 0,
1199: DM_POLYTOPE_TRIANGLE, 1, 0, 2, DM_POLYTOPE_TRIANGLE, 1, 1, 2, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 3, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 4, 0,
1200: DM_POLYTOPE_TRIANGLE, 1, 0, 3, DM_POLYTOPE_TRIANGLE, 1, 1, 3, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 2};
1201: static PetscInt ttriO[] = {0, 0, 0, 0,
1202: 0, 0, 0, 0,
1203: 0, 0, 0, 0,
1204: 0, 0, 0, -1, 0,
1205: 0, 0, 0, 0, -1,
1206: 0, 0, -1, 0, 0,
1207: 0, 0, 0, 0, 0};
1208: /* Add 1 edge and 4 tensor quads inside every tensor quad prism, making 4 new prisms. */
1209: static DMPolytopeType tquadT[] = {DM_POLYTOPE_POINT_PRISM_TENSOR, DM_POLYTOPE_SEG_PRISM_TENSOR, DM_POLYTOPE_QUAD_PRISM_TENSOR};
1210: static PetscInt tquadS[] = {1, 4, 4};
1211: static PetscInt tquadC[] = {DM_POLYTOPE_POINT, 1, 0, 0, DM_POLYTOPE_POINT, 1, 1, 0,
1212: DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 2, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 0, 0,
1213: DM_POLYTOPE_SEGMENT, 1, 0, 1, DM_POLYTOPE_SEGMENT, 1, 1, 1, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 3, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 0, 0,
1214: DM_POLYTOPE_SEGMENT, 1, 0, 2, DM_POLYTOPE_SEGMENT, 1, 1, 2, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 4, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 0, 0,
1215: DM_POLYTOPE_SEGMENT, 1, 0, 3, DM_POLYTOPE_SEGMENT, 1, 1, 3, DM_POLYTOPE_POINT_PRISM_TENSOR, 1, 5, 0, DM_POLYTOPE_POINT_PRISM_TENSOR, 0, 0,
1216: DM_POLYTOPE_QUADRILATERAL, 1, 0, 0, DM_POLYTOPE_QUADRILATERAL, 1, 1, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 2, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 3, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 5, 1,
1217: DM_POLYTOPE_QUADRILATERAL, 1, 0, 1, DM_POLYTOPE_QUADRILATERAL, 1, 1, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 2, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 3, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 0,
1218: DM_POLYTOPE_QUADRILATERAL, 1, 0, 2, DM_POLYTOPE_QUADRILATERAL, 1, 1, 2, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 3, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 4, 0, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 2,
1219: DM_POLYTOPE_QUADRILATERAL, 1, 0, 3, DM_POLYTOPE_QUADRILATERAL, 1, 1, 3, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 3, DM_POLYTOPE_SEG_PRISM_TENSOR, 0, 2, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 4, 1, DM_POLYTOPE_SEG_PRISM_TENSOR, 1, 5, 0};
1220: static PetscInt tquadO[] = {0, 0,
1221: 0, 0, 0, 0,
1222: 0, 0, 0, 0,
1223: 0, 0, 0, 0,
1224: 0, 0, 0, 0,
1225: 0, 0, 0, 0, -1, 0,
1226: 0, 0, 0, 0, 0, -1,
1227: 0, 0, -1, 0, 0, 0,
1228: 0, 0, 0, -1, 0, 0};
1229: /* Add 4 edges, 12 triangles, 1 quad, 4 tetrahedra, and 6 pyramids inside every pyramid. */
1230: static DMPolytopeType tpyrT[] = {DM_POLYTOPE_SEGMENT, DM_POLYTOPE_TRIANGLE, DM_POLYTOPE_QUADRILATERAL, DM_POLYTOPE_TETRAHEDRON, DM_POLYTOPE_PYRAMID};
1231: static PetscInt tpyrS[] = {4, 12, 1, 4, 6};
1232: static PetscInt tpyrC[] = {DM_POLYTOPE_POINT, 2, 1, 1, 0, DM_POLYTOPE_POINT, 1, 0, 0,
1233: DM_POLYTOPE_POINT, 2, 2, 1, 0, DM_POLYTOPE_POINT, 1, 0, 0,
1234: DM_POLYTOPE_POINT, 2, 3, 1, 0, DM_POLYTOPE_POINT, 1, 0, 0,
1235: DM_POLYTOPE_POINT, 2, 4, 1, 0, DM_POLYTOPE_POINT, 1, 0, 0,
1236: /* These four triangle face out of the bottom pyramid */
1237: DM_POLYTOPE_SEGMENT, 1, 1, 1, DM_POLYTOPE_SEGMENT, 0, 3, DM_POLYTOPE_SEGMENT, 0, 0,
1238: DM_POLYTOPE_SEGMENT, 1, 2, 1, DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 0, 1,
1239: DM_POLYTOPE_SEGMENT, 1, 3, 1, DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 0, 2,
1240: DM_POLYTOPE_SEGMENT, 1, 4, 1, DM_POLYTOPE_SEGMENT, 0, 2, DM_POLYTOPE_SEGMENT, 0, 3,
1241: /* These eight triangles face out of the corner pyramids */
1242: DM_POLYTOPE_SEGMENT, 1, 0, 3, DM_POLYTOPE_SEGMENT, 0, 3, DM_POLYTOPE_SEGMENT, 1, 1, 2,
1243: DM_POLYTOPE_SEGMENT, 1, 0, 2, DM_POLYTOPE_SEGMENT, 0, 0, DM_POLYTOPE_SEGMENT, 1, 2, 2,
1244: DM_POLYTOPE_SEGMENT, 1, 0, 1, DM_POLYTOPE_SEGMENT, 0, 1, DM_POLYTOPE_SEGMENT, 1, 3, 2,
1245: DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 0, 2, DM_POLYTOPE_SEGMENT, 1, 4, 2,
1246: DM_POLYTOPE_SEGMENT, 1, 0, 3, DM_POLYTOPE_SEGMENT, 1, 1, 0, DM_POLYTOPE_SEGMENT, 0, 0,
1247: DM_POLYTOPE_SEGMENT, 1, 0, 2, DM_POLYTOPE_SEGMENT, 1, 2, 0, DM_POLYTOPE_SEGMENT, 0, 1,
1248: DM_POLYTOPE_SEGMENT, 1, 0, 1, DM_POLYTOPE_SEGMENT, 1, 3, 0, DM_POLYTOPE_SEGMENT, 0, 2,
1249: DM_POLYTOPE_SEGMENT, 1, 0, 0, DM_POLYTOPE_SEGMENT, 1, 4, 0, DM_POLYTOPE_SEGMENT, 0, 3,
1250: /* This quad faces out of the bottom pyramid */
1251: DM_POLYTOPE_SEGMENT, 1, 1, 1, DM_POLYTOPE_SEGMENT, 1, 2, 1, DM_POLYTOPE_SEGMENT, 1, 3, 1, DM_POLYTOPE_SEGMENT, 1, 4, 1,
1252: /* The bottom face of each tet is on the triangular face */
1253: DM_POLYTOPE_TRIANGLE, 1, 1, 3, DM_POLYTOPE_TRIANGLE, 0, 8, DM_POLYTOPE_TRIANGLE, 0, 4, DM_POLYTOPE_TRIANGLE, 0, 0,
1254: DM_POLYTOPE_TRIANGLE, 1, 2, 3, DM_POLYTOPE_TRIANGLE, 0, 9, DM_POLYTOPE_TRIANGLE, 0, 5, DM_POLYTOPE_TRIANGLE, 0, 1,
1255: DM_POLYTOPE_TRIANGLE, 1, 3, 3, DM_POLYTOPE_TRIANGLE, 0, 10, DM_POLYTOPE_TRIANGLE, 0, 6, DM_POLYTOPE_TRIANGLE, 0, 2,
1256: DM_POLYTOPE_TRIANGLE, 1, 4, 3, DM_POLYTOPE_TRIANGLE, 0, 11, DM_POLYTOPE_TRIANGLE, 0, 7, DM_POLYTOPE_TRIANGLE, 0, 3,
1257: /* The front face of all pyramids is toward the front */
1258: DM_POLYTOPE_QUADRILATERAL, 1, 0, 0, DM_POLYTOPE_TRIANGLE, 1, 1, 0, DM_POLYTOPE_TRIANGLE, 0, 4, DM_POLYTOPE_TRIANGLE, 0, 11, DM_POLYTOPE_TRIANGLE, 1, 4, 1,
1259: DM_POLYTOPE_QUADRILATERAL, 1, 0, 3, DM_POLYTOPE_TRIANGLE, 1, 1, 1, DM_POLYTOPE_TRIANGLE, 1, 2, 0, DM_POLYTOPE_TRIANGLE, 0, 5, DM_POLYTOPE_TRIANGLE, 0, 8,
1260: DM_POLYTOPE_QUADRILATERAL, 1, 0, 2, DM_POLYTOPE_TRIANGLE, 0, 9, DM_POLYTOPE_TRIANGLE, 1, 2, 1, DM_POLYTOPE_TRIANGLE, 1, 3, 0, DM_POLYTOPE_TRIANGLE, 0, 6,
1261: DM_POLYTOPE_QUADRILATERAL, 1, 0, 1, DM_POLYTOPE_TRIANGLE, 0, 7, DM_POLYTOPE_TRIANGLE, 0, 10, DM_POLYTOPE_TRIANGLE, 1, 3, 1, DM_POLYTOPE_TRIANGLE, 1, 4, 0,
1262: DM_POLYTOPE_QUADRILATERAL, 0, 0, DM_POLYTOPE_TRIANGLE, 1, 1, 2, DM_POLYTOPE_TRIANGLE, 1, 2, 2, DM_POLYTOPE_TRIANGLE, 1, 3, 2, DM_POLYTOPE_TRIANGLE, 1, 4, 2,
1263: DM_POLYTOPE_QUADRILATERAL, 0, 0, DM_POLYTOPE_TRIANGLE, 0, 0, DM_POLYTOPE_TRIANGLE, 0, 3, DM_POLYTOPE_TRIANGLE, 0, 2, DM_POLYTOPE_TRIANGLE, 0, 1,
1264: };
1265: static PetscInt tpyrO[] = {0, 0,
1266: 0, 0,
1267: 0, 0,
1268: 0, 0,
1269: 0, 0, -1,
1270: 0, 0, -1,
1271: 0, 0, -1,
1272: 0, 0, -1,
1273: 0, -1, 0,
1274: 0, -1, 0,
1275: 0, -1, 0,
1276: 0, -1, 0,
1277: -1, 0, 0,
1278: -1, 0, 0,
1279: -1, 0, 0,
1280: -1, 0, 0,
1281: -1, -1, -1, -1,
1282: 0, -3, -2, -3,
1283: 0, -3, -2, -3,
1284: 0, -3, -2, -3,
1285: 0, -3, -2, -3,
1286: 0, 0, 0, 0, 0,
1287: 0, 0, 0, 0, 0,
1288: 0, 0, 0, 0, 0,
1289: 0, 0, 0, 0, 0,
1290: -2, 0, 0, 0, 0,
1291: 1, 0, 0, 0, 0};
1293: if (rt) *rt = 0;
1294: switch (source) {
1295: case DM_POLYTOPE_POINT: *Nt = 1; *target = vertexT; *size = vertexS; *cone = vertexC; *ornt = vertexO; break;
1296: case DM_POLYTOPE_SEGMENT: *Nt = 2; *target = segT; *size = segS; *cone = segC; *ornt = segO; break;
1297: case DM_POLYTOPE_POINT_PRISM_TENSOR: *Nt = 1; *target = tvertT; *size = tvertS; *cone = tvertC; *ornt = tvertO; break;
1298: case DM_POLYTOPE_TRIANGLE: *Nt = 2; *target = triT; *size = triS; *cone = triC; *ornt = triO; break;
1299: case DM_POLYTOPE_QUADRILATERAL: *Nt = 3; *target = quadT; *size = quadS; *cone = quadC; *ornt = quadO; break;
1300: case DM_POLYTOPE_SEG_PRISM_TENSOR: *Nt = 2; *target = tsegT; *size = tsegS; *cone = tsegC; *ornt = tsegO; break;
1301: case DM_POLYTOPE_TETRAHEDRON: *Nt = 3; *target = tetT; *size = tetS; *cone = tetC; *ornt = tetO; break;
1302: case DM_POLYTOPE_HEXAHEDRON: *Nt = 4; *target = hexT; *size = hexS; *cone = hexC; *ornt = hexO; break;
1303: case DM_POLYTOPE_TRI_PRISM: *Nt = 4; *target = tripT; *size = tripS; *cone = tripC; *ornt = tripO; break;
1304: case DM_POLYTOPE_TRI_PRISM_TENSOR: *Nt = 2; *target = ttriT; *size = ttriS; *cone = ttriC; *ornt = ttriO; break;
1305: case DM_POLYTOPE_QUAD_PRISM_TENSOR: *Nt = 3; *target = tquadT; *size = tquadS; *cone = tquadC; *ornt = tquadO; break;
1306: case DM_POLYTOPE_PYRAMID: *Nt = 5; *target = tpyrT; *size = tpyrS; *cone = tpyrC; *ornt = tpyrO; break;
1307: default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No refinement strategy for %s", DMPolytopeTypes[source]);
1308: }
1309: return 0;
1310: }
1312: static PetscErrorCode DMPlexTransformInitialize_Regular(DMPlexTransform tr)
1313: {
1314: tr->ops->view = DMPlexTransformView_Regular;
1315: tr->ops->setup = DMPlexTransformSetUp_Regular;
1316: tr->ops->destroy = DMPlexTransformDestroy_Regular;
1317: tr->ops->celltransform = DMPlexTransformCellRefine_Regular;
1318: tr->ops->getsubcellorientation = DMPlexTransformGetSubcellOrientation_Regular;
1319: tr->ops->mapcoordinates = DMPlexTransformMapCoordinatesBarycenter_Internal;
1320: return 0;
1321: }
1323: PETSC_EXTERN PetscErrorCode DMPlexTransformCreate_Regular(DMPlexTransform tr)
1324: {
1325: DMPlexRefine_Regular *f;
1328: PetscNewLog(tr, &f);
1329: tr->data = f;
1331: DMPlexTransformInitialize_Regular(tr);
1332: return 0;
1333: }