xref: /petsc/src/dm/impls/plex/tests/ex37.c (revision 5b72a2ce2e245e80f3c4d5fd294a594e197f5e74)
1c4762a1bSJed Brown static const char help[] = "Test of EGADSLite CAD functionality";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <egads.h>
6c4762a1bSJed Brown #include <petsc.h>
7c4762a1bSJed Brown 
8c4762a1bSJed Brown typedef struct {
9c4762a1bSJed Brown   char filename[PETSC_MAX_PATH_LEN];
10c4762a1bSJed Brown } AppCtx;
11c4762a1bSJed Brown 
12c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
13c4762a1bSJed Brown {
14c4762a1bSJed Brown   PetscErrorCode ierr;
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   PetscFunctionBeginUser;
17c4762a1bSJed Brown   options->filename[0] = '\0';
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "EGADSPlex Problem Options", "EGADSLite");CHKERRQ(ierr);
20c4762a1bSJed Brown   ierr = PetscOptionsString("-filename", "The EGADSLite file", "ex9.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
21c4762a1bSJed Brown   ierr = PetscOptionsEnd();
22c4762a1bSJed Brown   PetscFunctionReturn(0);
23c4762a1bSJed Brown }
24c4762a1bSJed Brown 
25c4762a1bSJed Brown int main(int argc, char *argv[])
26c4762a1bSJed Brown {
27c4762a1bSJed Brown   DMLabel        bodyLabel, faceLabel, edgeLabel;
28c4762a1bSJed Brown   PetscInt       cStart, cEnd, c;
29c4762a1bSJed Brown   /* EGADSLite variables */
30c4762a1bSJed Brown   ego            context, model, geom, *bodies, *objs, *nobjs, *mobjs, *lobjs;
31c4762a1bSJed Brown   int            oclass, mtype, nbodies, *senses;
32c4762a1bSJed Brown   int            b;
33c4762a1bSJed Brown   /* PETSc variables */
34c4762a1bSJed Brown   DM             dm;
35c4762a1bSJed Brown   PetscInt       dim = -1, cdim = -1, numCorners = 0, numVertices = 0, numCells = 0;
36c4762a1bSJed Brown   PetscInt      *cells  = NULL;
37c4762a1bSJed Brown   PetscReal     *coords = NULL;
38c4762a1bSJed Brown   MPI_Comm       comm;
39c4762a1bSJed Brown   PetscMPIInt    rank;
40c4762a1bSJed Brown   AppCtx         ctx;
41c4762a1bSJed Brown   PetscErrorCode ierr;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
44c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
45c4762a1bSJed Brown   ierr = ProcessOptions(comm, &ctx);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
47c4762a1bSJed Brown   if (!rank) {
48c4762a1bSJed Brown     /* Open EGADs file and load EGADs model data */
49c4762a1bSJed Brown     ierr = EG_open(&context);CHKERRQ(ierr);
50c4762a1bSJed Brown     ierr = EG_loadModel(context, 0, ctx.filename, &model);CHKERRQ(ierr);
51c4762a1bSJed Brown 
52c4762a1bSJed Brown     /* test bodyTopo functions */
53c4762a1bSJed Brown     ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses);CHKERRQ(ierr);
54c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, " Number of BODIES (nbodies): %d \n", nbodies);CHKERRQ(ierr);
55c4762a1bSJed Brown 
56c4762a1bSJed Brown     for (b = 0; b < nbodies; ++b) {
57c4762a1bSJed Brown       ego body = bodies[b];
58c4762a1bSJed Brown       int id, Nsh, Nf, Nl, l, Ne, e, Nv, v;
59c4762a1bSJed Brown 
60c4762a1bSJed Brown       /* Output Basic Model Topology */
61c4762a1bSJed Brown       ierr = EG_getBodyTopos(body, NULL, SHELL, &Nsh, &objs);CHKERRQ(ierr);
62c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of SHELLS: %d \n", Nsh);CHKERRQ(ierr);
63c4762a1bSJed Brown 
64c4762a1bSJed Brown       ierr = EG_getBodyTopos(body, NULL, FACE,  &Nf, &objs);CHKERRQ(ierr);
65c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of FACES: %d \n", Nf);CHKERRQ(ierr);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown       ierr = EG_getBodyTopos(body, NULL, LOOP,  &Nl, &lobjs);CHKERRQ(ierr);
68c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of LOOPS: %d \n", Nl);CHKERRQ(ierr);
69c4762a1bSJed Brown 
70c4762a1bSJed Brown       ierr = EG_getBodyTopos(body, NULL, EDGE,  &Ne, &objs);CHKERRQ(ierr);
71c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of EDGES: %d \n", Ne);CHKERRQ(ierr);
72c4762a1bSJed Brown 
73c4762a1bSJed Brown       ierr = EG_getBodyTopos(body, NULL, NODE,  &Nv, &objs);CHKERRQ(ierr);
74c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of NODES: %d \n", Nv);CHKERRQ(ierr);
75c4762a1bSJed Brown 
76c4762a1bSJed Brown       for (l = 0; l < Nl; ++l) {
77c4762a1bSJed Brown         ego loop = lobjs[l];
78c4762a1bSJed Brown 
79c4762a1bSJed Brown         id   = EG_indexBodyTopo(body, loop);
80c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "          LOOP ID: %d\n", id);CHKERRQ(ierr);
81c4762a1bSJed Brown 
82c4762a1bSJed Brown         /* Get EDGE info which associated with the current LOOP */
83c4762a1bSJed Brown         ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
84c4762a1bSJed Brown 
85c4762a1bSJed Brown         for (e = 0; e < Ne; ++e) {
86c4762a1bSJed Brown           ego edge = objs[e];
87c4762a1bSJed Brown 
88c4762a1bSJed Brown           id = EG_indexBodyTopo(body, edge);CHKERRQ(ierr);
89c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_SELF, "            EDGE ID: %d\n", id);CHKERRQ(ierr);
90c4762a1bSJed Brown 
91c4762a1bSJed Brown           double range[4] = {0., 0., 0., 0.};
92c4762a1bSJed Brown           double point[3] = {0., 0., 0.};
93c4762a1bSJed Brown           int    peri;
94c4762a1bSJed Brown 
95c4762a1bSJed Brown           ierr = EG_getRange(objs[e], range, &peri);
96c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_SELF, " Range = %lf, %lf, %lf, %lf \n", range[0], range[1], range[2], range[3]);
97c4762a1bSJed Brown 
98c4762a1bSJed Brown           /* Get NODE info which associated with the current EDGE */
99c4762a1bSJed Brown           ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);CHKERRQ(ierr);
100c4762a1bSJed Brown 
101c4762a1bSJed Brown           for (v = 0; v < Nv; ++v) {
102c4762a1bSJed Brown             ego    vertex = nobjs[v];
103c4762a1bSJed Brown             double limits[4];
104c4762a1bSJed Brown             int    dummy;
105c4762a1bSJed Brown 
106c4762a1bSJed Brown             ierr = EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses);CHKERRQ(ierr);
107c4762a1bSJed Brown             id   = EG_indexBodyTopo(body, vertex);
108c4762a1bSJed Brown             ierr = PetscPrintf(PETSC_COMM_SELF, "              NODE ID: %d \n", id);CHKERRQ(ierr);
109c4762a1bSJed Brown             ierr = PetscPrintf(PETSC_COMM_SELF, "                 (x, y, z) = (%lf, %lf, %lf) \n", limits[0], limits[1], limits[2]);
110c4762a1bSJed Brown 
111c4762a1bSJed Brown             point[0] = point[0] + limits[0];
112c4762a1bSJed Brown             point[1] = point[1] + limits[1];
113c4762a1bSJed Brown             point[2] = point[2] + limits[2];
114c4762a1bSJed Brown           }
115c4762a1bSJed Brown         }
116c4762a1bSJed Brown       }
117c4762a1bSJed Brown     }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown     /* ---------------------------------------------------------------------------------------------------
120c4762a1bSJed Brown     Generate Petsc Plex
121c4762a1bSJed Brown       Get all Nodes in model, record coordinates in a correctly formatted array
122c4762a1bSJed Brown       Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array */
123c4762a1bSJed Brown 
124c4762a1bSJed Brown     /* Calculate cell and vertex sizes */
125c4762a1bSJed Brown     ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses);CHKERRQ(ierr);
126c4762a1bSJed Brown     numCells    = 0;
127c4762a1bSJed Brown     numVertices = 0;
128c4762a1bSJed Brown     for (b = 0; b < nbodies; ++b) {
129c4762a1bSJed Brown       ego body = bodies[b];
130c4762a1bSJed Brown       int id, Nl, l, Nv, v;
131c4762a1bSJed Brown 
132c4762a1bSJed Brown       ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
133c4762a1bSJed Brown       ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs);CHKERRQ(ierr);
134c4762a1bSJed Brown       for (l = 0; l < Nl; ++l) {
135c4762a1bSJed Brown         ego loop = lobjs[l];
136c4762a1bSJed Brown 
137c4762a1bSJed Brown         id = EG_indexBodyTopo(body, loop);
138c4762a1bSJed Brown         /* TODO: Instead of assuming contiguous ids, we could use a hash table */
139c4762a1bSJed Brown         numCells = PetscMax(id, numCells);
140c4762a1bSJed Brown       }
141c4762a1bSJed Brown       for (v = 0; v < Nv; ++v) {
142c4762a1bSJed Brown         ego vertex = nobjs[v];
143c4762a1bSJed Brown 
144c4762a1bSJed Brown         id = EG_indexBodyTopo(body, vertex);
145c4762a1bSJed Brown         /* TODO: Instead of assuming contiguous ids, we could use a hash table */
146c4762a1bSJed Brown         numVertices = PetscMax(id, numVertices);
147c4762a1bSJed Brown       }
148c4762a1bSJed Brown     }
149c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, "\nPLEX Input Array Checkouts\n");CHKERRQ(ierr);
150c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, " Total Number of Unique Cells    = %d \n", numCells);CHKERRQ(ierr);
151c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF, " Total Number of Unique Vertices = %d \n", numVertices);CHKERRQ(ierr);
152c4762a1bSJed Brown 
153c4762a1bSJed Brown     dim        = 2; /* Assume 3D Models :: Need to update to handle 2D Models in the future */
154c4762a1bSJed Brown     cdim       = 3; /* Assume 3D Models :: Need to update to handle 2D Models in the future */
155c4762a1bSJed Brown     numCorners = 3; /* TODO Check number of cell corners from EGADSLite */
156c4762a1bSJed Brown     ierr = PetscMalloc2(numVertices*cdim, &coords, numCells*numCorners, &cells);CHKERRQ(ierr);
157c4762a1bSJed Brown 
158c4762a1bSJed Brown     /* Get vertex coordinates */
159c4762a1bSJed Brown     for (b = 0; b < nbodies; ++b) {
160c4762a1bSJed Brown       ego body = bodies[b];
161c4762a1bSJed Brown       int id, Nv, v;
162c4762a1bSJed Brown 
163c4762a1bSJed Brown       ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs);CHKERRQ(ierr);
164c4762a1bSJed Brown       for (v = 0; v < Nv; ++v) {
165c4762a1bSJed Brown         ego    vertex = nobjs[v];
166c4762a1bSJed Brown         double limits[4];
167c4762a1bSJed Brown         int    dummy;
168c4762a1bSJed Brown 
169c4762a1bSJed Brown         ierr = EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses);CHKERRQ(ierr);
170c4762a1bSJed Brown         id   = EG_indexBodyTopo(body, vertex);CHKERRQ(ierr);
171c4762a1bSJed Brown         coords[(id-1)*cdim+0] = limits[0];
172c4762a1bSJed Brown         coords[(id-1)*cdim+1] = limits[1];
173c4762a1bSJed Brown         coords[(id-1)*cdim+2] = limits[2];
174c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "    Node ID = %d \n", id);
175c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "      (x,y,z) = (%lf, %lf, %lf) \n \n", coords[(id-1)*cdim+0], coords[(id-1)*cdim+1],coords[(id-1)*cdim+2]);
176c4762a1bSJed Brown       }
177c4762a1bSJed Brown     }
178c4762a1bSJed Brown 
179c4762a1bSJed Brown     /* Get cell vertices by traversing loops */
180c4762a1bSJed Brown     for (b = 0; b < nbodies; ++b) {
181c4762a1bSJed Brown       ego body = bodies[b];
182c4762a1bSJed Brown       int id, Nl, l;
183c4762a1bSJed Brown 
184c4762a1bSJed Brown       ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
185c4762a1bSJed Brown       for (l = 0; l < Nl; ++l) {
186c4762a1bSJed Brown         ego loop = lobjs[l];
187c4762a1bSJed Brown         int lid, Ne, e, nc = 0, c;
188c4762a1bSJed Brown 
189c4762a1bSJed Brown         lid  = EG_indexBodyTopo(body, loop);CHKERRQ(ierr);
190c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "    LOOP ID: %d \n", lid);CHKERRQ(ierr);
191c4762a1bSJed Brown         ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
192c4762a1bSJed Brown 
193c4762a1bSJed Brown         for (e = 0; e < Ne; ++e) {
194c4762a1bSJed Brown           ego edge = objs[e];
195c4762a1bSJed Brown           int Nv, v;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown           id   = EG_indexBodyTopo(body, edge);
198c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_SELF, "      EDGE ID: %d \n", id);CHKERRQ(ierr);
199c4762a1bSJed Brown           if (mtype == DEGENERATE) {ierr = PetscPrintf(PETSC_COMM_SELF, "        EGDE %d is DEGENERATE \n", id);CHKERRQ(ierr);}
200c4762a1bSJed Brown           ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);
201c4762a1bSJed Brown 
202c4762a1bSJed Brown           /* Add unique vertices to cells, this handles mtype == DEGENERATE fine */
203c4762a1bSJed Brown           for (v = 0; v < Nv; ++v) {
204c4762a1bSJed Brown             ego vertex = nobjs[v];
205c4762a1bSJed Brown 
206c4762a1bSJed Brown             id = EG_indexBodyTopo(body, vertex);
207c4762a1bSJed Brown             for (c = 0; c < nc; ++c) if (cells[(lid-1)*numCorners+c] == id-1) break;
208c4762a1bSJed Brown             if (c == nc) cells[(lid-1)*numCorners+nc++] = id-1;
209c4762a1bSJed Brown           }
210c4762a1bSJed Brown         }
211c4762a1bSJed Brown         if (nc != numCorners) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid number of cell corners %D, should be %D", nc, numCorners);
212c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, "      LOOP Corner NODEs (");
213c4762a1bSJed Brown         for (c = 0; c < numCorners; ++c) {
214c4762a1bSJed Brown           if (c > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, ", ");}
215c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_SELF, "%D", cells[(lid-1)*numCorners+c]);
216c4762a1bSJed Brown         }
217c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF, ")\n");
218c4762a1bSJed Brown       }
219c4762a1bSJed Brown     }
220c4762a1bSJed Brown   }
221c4762a1bSJed Brown   ierr = DMPlexCreateFromCellList(PETSC_COMM_WORLD, dim, numCells, numVertices, numCorners, PETSC_TRUE, cells, cdim, coords, &dm);CHKERRQ(ierr);
222c4762a1bSJed Brown   ierr = PetscFree2(coords, cells);CHKERRQ(ierr);
223c4762a1bSJed Brown   {
224c4762a1bSJed Brown     PetscContainer modelObj;
225c4762a1bSJed Brown 
226c4762a1bSJed Brown     ierr = PetscContainerCreate(PETSC_COMM_SELF, &modelObj);CHKERRQ(ierr);
227c4762a1bSJed Brown     ierr = PetscContainerSetPointer(modelObj, model);CHKERRQ(ierr);
228c4762a1bSJed Brown     ierr = PetscObjectCompose((PetscObject) dm, "EGADS Model", (PetscObject) modelObj);CHKERRQ(ierr);
229c4762a1bSJed Brown     ierr = PetscContainerDestroy(&modelObj);CHKERRQ(ierr);
230c4762a1bSJed Brown   }
231c4762a1bSJed Brown   ierr = DMCreateLabel(dm, "EGADS Body ID");CHKERRQ(ierr);
232c4762a1bSJed Brown   ierr = DMGetLabel(dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
233c4762a1bSJed Brown   ierr = DMCreateLabel(dm, "EGADS Face ID");CHKERRQ(ierr);
234c4762a1bSJed Brown   ierr = DMGetLabel(dm, "EGADS Face ID", &faceLabel);CHKERRQ(ierr);
235c4762a1bSJed Brown   ierr = DMCreateLabel(dm, "EGADS Edge ID");CHKERRQ(ierr);
236c4762a1bSJed Brown   ierr = DMGetLabel(dm, "EGADS Edge ID", &edgeLabel);CHKERRQ(ierr);
237c4762a1bSJed Brown   for (b = 0; b < nbodies; ++b) {
238c4762a1bSJed Brown     ego body = bodies[b];
239c4762a1bSJed Brown     int id, Nl, l;
240c4762a1bSJed Brown 
241c4762a1bSJed Brown     ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
242c4762a1bSJed Brown     for (l = 0; l < Nl; ++l) {
243c4762a1bSJed Brown       ego loop = lobjs[l];
244c4762a1bSJed Brown       int lid, cell, Ne, e;
245c4762a1bSJed Brown 
246c4762a1bSJed Brown       lid  = EG_indexBodyTopo(body, loop);CHKERRQ(ierr);
247c4762a1bSJed Brown       cell = lid-1;
248c4762a1bSJed Brown       ierr = DMLabelSetValue(bodyLabel, cell, b);CHKERRQ(ierr);
249c4762a1bSJed Brown       {
250c4762a1bSJed Brown         ego *fobjs;
251c4762a1bSJed Brown         int  Nf, fid;
252c4762a1bSJed Brown 
253c4762a1bSJed Brown         ierr = EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs);CHKERRQ(ierr);
254c4762a1bSJed Brown         fid  = EG_indexBodyTopo(body, fobjs[0]);CHKERRQ(ierr);
255c4762a1bSJed Brown         ierr = DMLabelSetValue(faceLabel, cell, fid);CHKERRQ(ierr);
256c4762a1bSJed Brown       }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown       ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
259c4762a1bSJed Brown       for (e = 0; e < Ne; ++e) {
260c4762a1bSJed Brown         ego             edge = objs[e];
261c4762a1bSJed Brown         int             eid, Nv, v;
262c4762a1bSJed Brown         PetscInt        support[2], numEdges;
263c4762a1bSJed Brown         const PetscInt *edges;
264c4762a1bSJed Brown 
265c4762a1bSJed Brown         eid  = EG_indexBodyTopo(body, edge);
266c4762a1bSJed Brown         ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);
267c4762a1bSJed Brown         if (Nv > 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Edge %d has %d vertices > 2", eid, Nv);
268c4762a1bSJed Brown         for (v = 0; v < Nv; ++v) {
269c4762a1bSJed Brown           ego vertex = nobjs[v];
270c4762a1bSJed Brown 
271c4762a1bSJed Brown           id   = EG_indexBodyTopo(body, vertex);
272c4762a1bSJed Brown           ierr = DMLabelSetValue(edgeLabel, numCells + id-1, eid);CHKERRQ(ierr);
273c4762a1bSJed Brown           support[v] = numCells + id-1;
274c4762a1bSJed Brown         }
275c4762a1bSJed Brown         if (Nv == 2) {
276c4762a1bSJed Brown           ierr = DMPlexGetJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr);
277c4762a1bSJed Brown           if (numEdges != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "2 vertices should only bound 1 edge, not %D", numEdges);
278c4762a1bSJed Brown           ierr = DMLabelSetValue(edgeLabel, edges[0], eid);CHKERRQ(ierr);
279c4762a1bSJed Brown           ierr = DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr);
280c4762a1bSJed Brown         }
281c4762a1bSJed Brown       }
282c4762a1bSJed Brown     }
283c4762a1bSJed Brown   }
284c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
285c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
286c4762a1bSJed Brown     PetscInt *closure = NULL;
287c4762a1bSJed Brown     PetscInt  clSize, cl, bval, fval;
288c4762a1bSJed Brown 
289c4762a1bSJed Brown     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
290c4762a1bSJed Brown     ierr = DMLabelGetValue(bodyLabel, c, &bval);CHKERRQ(ierr);
291c4762a1bSJed Brown     ierr = DMLabelGetValue(faceLabel, c, &fval);CHKERRQ(ierr);
292c4762a1bSJed Brown     for (cl = 0; cl < clSize*2; cl += 2) {
293c4762a1bSJed Brown       ierr = DMLabelSetValue(bodyLabel, closure[cl], bval);CHKERRQ(ierr);
294c4762a1bSJed Brown       ierr = DMLabelSetValue(faceLabel, closure[cl], fval);CHKERRQ(ierr);
295c4762a1bSJed Brown     }
296c4762a1bSJed Brown     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
297c4762a1bSJed Brown   }
298c4762a1bSJed Brown   ierr = DMLabelView(bodyLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
299c4762a1bSJed Brown   ierr = DMLabelView(faceLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
300c4762a1bSJed Brown   ierr = DMLabelView(edgeLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
301c4762a1bSJed Brown   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr);
304c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   /* Close EGADSlite file */
307c4762a1bSJed Brown   ierr = EG_close(context);CHKERRQ(ierr);
308c4762a1bSJed Brown   ierr = PetscFinalize();
309c4762a1bSJed Brown   return ierr;
310c4762a1bSJed Brown }
311c4762a1bSJed Brown 
312c4762a1bSJed Brown /*TEST
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   build:
315c4762a1bSJed Brown     requires: egads
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   test:
318c4762a1bSJed Brown     suffix: sphere_0
319*5b72a2ceSPierre Jolivet     filter: sed "s/DM_[0-9a-zA-Z]*_0/DM__0/g"
320c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/unit_sphere.egadslite -dm_view ::ascii_info_detail
321c4762a1bSJed Brown 
322c4762a1bSJed Brown TEST*/
323