xref: /petsc/src/dm/impls/plex/plexegads.c (revision c03d1177128e9544fa7c6deef245a4c5fa59c988)
1a8ededdfSMatthew G. Knepley #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
27bee2925SMatthew Knepley #include <petsc/private/hashmapi.h>
3a8ededdfSMatthew G. Knepley 
4a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS
5a8ededdfSMatthew G. Knepley #include <egads.h>
6a8ededdfSMatthew G. Knepley #endif
7a8ededdfSMatthew G. Knepley 
8a8ededdfSMatthew G. Knepley /* We need to understand how to natively parse STEP files. There seems to be only one open source implementation of
9a8ededdfSMatthew G. Knepley    the STEP parser contained in the OpenCASCADE package. It is enough to make a strong man weep:
10a8ededdfSMatthew G. Knepley 
11a8ededdfSMatthew G. Knepley      https://github.com/tpaviot/oce/tree/master/src/STEPControl
12a8ededdfSMatthew G. Knepley 
13a8ededdfSMatthew G. Knepley    The STEP, and inner EXPRESS, formats are ISO standards, so they are documented
14a8ededdfSMatthew G. Knepley 
15a8ededdfSMatthew G. Knepley      https://stackoverflow.com/questions/26774037/documentation-or-specification-for-step-and-stp-files
16a8ededdfSMatthew G. Knepley      http://stepmod.sourceforge.net/express_model_spec/
17a8ededdfSMatthew G. Knepley 
18a8ededdfSMatthew G. Knepley    but again it seems that there has been a deliberate effort at obfuscation, probably to raise the bar for entrants.
19a8ededdfSMatthew G. Knepley */
20a8ededdfSMatthew G. Knepley 
21a8ededdfSMatthew G. Knepley 
22a8ededdfSMatthew G. Knepley /*@
23a8ededdfSMatthew G. Knepley   DMPlexSnapToGeomModel - Given a coordinate point 'mcoords' on the mesh point 'p', return the closest coordinate point 'gcoords' on the geometry model associated with that point.
24a8ededdfSMatthew G. Knepley 
25a8ededdfSMatthew G. Knepley   Not collective
26a8ededdfSMatthew G. Knepley 
27a8ededdfSMatthew G. Knepley   Input Parameters:
28a8ededdfSMatthew G. Knepley + dm      - The DMPlex object
29a8ededdfSMatthew G. Knepley . p       - The mesh point
30a8ededdfSMatthew G. Knepley - mcoords - A coordinate point lying on the mesh point
31a8ededdfSMatthew G. Knepley 
32a8ededdfSMatthew G. Knepley   Output Parameter:
33a8ededdfSMatthew G. Knepley . gcoords - The closest coordinate point on the geometry model associated with 'p' to the given point
34a8ededdfSMatthew G. Knepley 
35a8ededdfSMatthew G. Knepley   Note: Returns the original coordinates if no geometry model is found. Right now the only supported geometry model is EGADS.
36a8ededdfSMatthew G. Knepley 
37a8ededdfSMatthew G. Knepley   Level: intermediate
38a8ededdfSMatthew G. Knepley 
39a8ededdfSMatthew G. Knepley .seealso: DMRefine(), DMPlexCreate(), DMPlexSetRefinementUniform()
40a8ededdfSMatthew G. Knepley @*/
41a8ededdfSMatthew G. Knepley PetscErrorCode DMPlexSnapToGeomModel(DM dm, PetscInt p, const PetscScalar mcoords[], PetscScalar gcoords[])
42a8ededdfSMatthew G. Knepley {
43a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS
44f0fcf0adSMatthew G. Knepley   DM             cdm;
45a8ededdfSMatthew G. Knepley   DMLabel        bodyLabel, faceLabel, edgeLabel;
46a8ededdfSMatthew G. Knepley   PetscContainer modelObj;
47a8ededdfSMatthew G. Knepley   PetscInt       bodyID, faceID, edgeID;
48a8ededdfSMatthew G. Knepley   ego           *bodies;
49f0fcf0adSMatthew G. Knepley   ego            model, geom, body, obj;
50f0fcf0adSMatthew G. Knepley   /* result has to hold derviatives, along with the value */
51f0fcf0adSMatthew G. Knepley   double         params[3], result[18], paramsV[16*3], resultV[16*3];
52a8ededdfSMatthew G. Knepley   int            Nb, oclass, mtype, *senses;
53f0fcf0adSMatthew G. Knepley   Vec            coordinatesLocal;
54f0fcf0adSMatthew G. Knepley   PetscScalar   *coords = NULL;
55f0fcf0adSMatthew G. Knepley   PetscInt       Nv, v, Np = 0, pm;
56a8ededdfSMatthew G. Knepley #endif
57f0fcf0adSMatthew G. Knepley   PetscInt       dE, d;
58a8ededdfSMatthew G. Knepley   PetscErrorCode ierr;
59a8ededdfSMatthew G. Knepley 
60a8ededdfSMatthew G. Knepley   PetscFunctionBegin;
61f0fcf0adSMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr);
62a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS
63a8ededdfSMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Body ID",   &bodyLabel);CHKERRQ(ierr);
64a8ededdfSMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Face ID",   &faceLabel);CHKERRQ(ierr);
65a8ededdfSMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Edge ID",   &edgeLabel);CHKERRQ(ierr);
66a8ededdfSMatthew G. Knepley   if (!bodyLabel || !faceLabel || !edgeLabel) {
67f0fcf0adSMatthew G. Knepley     for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
68a8ededdfSMatthew G. Knepley     PetscFunctionReturn(0);
69a8ededdfSMatthew G. Knepley   }
70f0fcf0adSMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
71f0fcf0adSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinatesLocal);CHKERRQ(ierr);
72a8ededdfSMatthew G. Knepley   ierr = DMLabelGetValue(bodyLabel,   p, &bodyID);CHKERRQ(ierr);
73a8ededdfSMatthew G. Knepley   ierr = DMLabelGetValue(faceLabel,   p, &faceID);CHKERRQ(ierr);
74a8ededdfSMatthew G. Knepley   ierr = DMLabelGetValue(edgeLabel,   p, &edgeID);CHKERRQ(ierr);
75a8ededdfSMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
76a8ededdfSMatthew G. Knepley   ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr);
77a8ededdfSMatthew G. Knepley   ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
78f0fcf0adSMatthew G. Knepley   if (bodyID < 0 || bodyID >= Nb) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %D is not in [0, %d)", bodyID, Nb);
79a8ededdfSMatthew G. Knepley   body = bodies[bodyID];
80f0fcf0adSMatthew G. Knepley 
81f0fcf0adSMatthew G. Knepley   if (edgeID >= 0)      {ierr = EG_objectBodyTopo(body, EDGE, edgeID, &obj);CHKERRQ(ierr); Np = 1;}
82f0fcf0adSMatthew G. Knepley   else if (faceID >= 0) {ierr = EG_objectBodyTopo(body, FACE, faceID, &obj);CHKERRQ(ierr); Np = 2;}
83f0fcf0adSMatthew G. Knepley   else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Point %D is not in edge or face label for EGADS", p);
84f0fcf0adSMatthew G. Knepley   /* Calculate parameters (t or u,v) for vertices */
85f0fcf0adSMatthew G. Knepley   ierr = DMPlexVecGetClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords);CHKERRQ(ierr);
86f0fcf0adSMatthew G. Knepley   Nv  /= dE;
87f0fcf0adSMatthew G. Knepley   if (Nv == 1) {
88f0fcf0adSMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords);CHKERRQ(ierr);
89f0fcf0adSMatthew G. Knepley     for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
90f0fcf0adSMatthew G. Knepley     PetscFunctionReturn(0);
91a8ededdfSMatthew G. Knepley   }
92f0fcf0adSMatthew G. Knepley   if (Nv > 16) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %D coordinates associated to point %D", Nv, p);
93f0fcf0adSMatthew G. Knepley   for (v = 0; v < Nv; ++v) {ierr = EG_invEvaluate(obj, &coords[v*dE], &paramsV[v*3], &resultV[v*3]);CHKERRQ(ierr);}
94f0fcf0adSMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords);CHKERRQ(ierr);
95f0fcf0adSMatthew G. Knepley   /* Calculate parameters (t or u,v) for new vertex at edge midpoint */
96f0fcf0adSMatthew G. Knepley   for (pm = 0; pm < Np; ++pm) {
97f0fcf0adSMatthew G. Knepley     params[pm] = 0.;
98f0fcf0adSMatthew G. Knepley     for (v = 0; v < Nv; ++v) {params[pm] += paramsV[v*3+pm];}
99f0fcf0adSMatthew G. Knepley     params[pm] /= Nv;
100f0fcf0adSMatthew G. Knepley   }
101f0fcf0adSMatthew G. Knepley   /* TODO Check
102f0fcf0adSMatthew G. Knepley     double range[4]; // [umin, umax, vmin, vmax]
103f0fcf0adSMatthew G. Knepley     int    peri;
104f0fcf0adSMatthew G. Knepley     ierr = EG_getRange(face, range, &peri);CHKERRQ(ierr);
105f0fcf0adSMatthew G. Knepley     if ((paramsNew[0] < range[0]) || (paramsNew[0] > range[1]) || (paramsNew[1] < range[2]) || (paramsNew[1] > range[3])) SETERRQ();
106f0fcf0adSMatthew G. Knepley   */
107f0fcf0adSMatthew G. Knepley   /* Put coordinates for new vertex in result[] */
108f0fcf0adSMatthew G. Knepley   ierr = EG_evaluate(obj, params, result);CHKERRQ(ierr);
109f0fcf0adSMatthew G. Knepley   for (d = 0; d < dE; ++d) gcoords[d] = result[d];
110a8ededdfSMatthew G. Knepley #else
111f0fcf0adSMatthew G. Knepley   for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
112a8ededdfSMatthew G. Knepley #endif
113a8ededdfSMatthew G. Knepley   PetscFunctionReturn(0);
114a8ededdfSMatthew G. Knepley }
1157bee2925SMatthew Knepley 
1167bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS)
1177bee2925SMatthew Knepley static PetscErrorCode DMPlexEGADSPrintModel(ego model)
1187bee2925SMatthew Knepley {
1197bee2925SMatthew Knepley   ego            geom, *bodies, *objs, *nobjs, *mobjs, *lobjs;
1207bee2925SMatthew Knepley   int            oclass, mtype, *senses;
1217bee2925SMatthew Knepley   int            Nb, b;
1227bee2925SMatthew Knepley   PetscErrorCode ierr;
1237bee2925SMatthew Knepley 
1247bee2925SMatthew Knepley   PetscFunctionBeginUser;
1257bee2925SMatthew Knepley   /* test bodyTopo functions */
1267bee2925SMatthew Knepley   ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
1277bee2925SMatthew Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, " Number of BODIES (nbodies): %d \n", Nb);CHKERRQ(ierr);
1287bee2925SMatthew Knepley 
1297bee2925SMatthew Knepley   for (b = 0; b < Nb; ++b) {
1307bee2925SMatthew Knepley     ego body = bodies[b];
1317bee2925SMatthew Knepley     int id, Nsh, Nf, Nl, l, Ne, e, Nv, v;
1327bee2925SMatthew Knepley 
1337bee2925SMatthew Knepley     /* Output Basic Model Topology */
1347bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, SHELL, &Nsh, &objs);CHKERRQ(ierr);
1357bee2925SMatthew Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of SHELLS: %d \n", Nsh);CHKERRQ(ierr);
1367bee2925SMatthew Knepley     EG_free(objs);
1377bee2925SMatthew Knepley 
1387bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, FACE,  &Nf, &objs);CHKERRQ(ierr);
1397bee2925SMatthew Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of FACES: %d \n", Nf);CHKERRQ(ierr);
1407bee2925SMatthew Knepley     EG_free(objs);
1417bee2925SMatthew Knepley 
1427bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, LOOP,  &Nl, &lobjs);CHKERRQ(ierr);
1437bee2925SMatthew Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of LOOPS: %d \n", Nl);CHKERRQ(ierr);
1447bee2925SMatthew Knepley 
1457bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, EDGE,  &Ne, &objs);CHKERRQ(ierr);
1467bee2925SMatthew Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of EDGES: %d \n", Ne);CHKERRQ(ierr);
1477bee2925SMatthew Knepley     EG_free(objs);
1487bee2925SMatthew Knepley 
1497bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, NODE,  &Nv, &objs);CHKERRQ(ierr);
1507bee2925SMatthew Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of NODES: %d \n", Nv);CHKERRQ(ierr);
1517bee2925SMatthew Knepley     EG_free(objs);
1527bee2925SMatthew Knepley 
1537bee2925SMatthew Knepley     for (l = 0; l < Nl; ++l) {
1547bee2925SMatthew Knepley       ego loop = lobjs[l];
1557bee2925SMatthew Knepley 
1567bee2925SMatthew Knepley       id   = EG_indexBodyTopo(body, loop);
1577bee2925SMatthew Knepley       ierr = PetscPrintf(PETSC_COMM_SELF, "          LOOP ID: %d\n", id);CHKERRQ(ierr);
1587bee2925SMatthew Knepley 
1597bee2925SMatthew Knepley       /* Get EDGE info which associated with the current LOOP */
1607bee2925SMatthew Knepley       ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
1617bee2925SMatthew Knepley 
1627bee2925SMatthew Knepley       for (e = 0; e < Ne; ++e) {
1637bee2925SMatthew Knepley         ego    edge      = objs[e];
1647bee2925SMatthew Knepley         double range[4]  = {0., 0., 0., 0.};
1657bee2925SMatthew Knepley         double point[3]  = {0., 0., 0.};
166266cfabeSMatthew G. Knepley         double params[3] = {0., 0., 0.};
167266cfabeSMatthew G. Knepley         double result[18];
1687bee2925SMatthew Knepley         int    peri;
1697bee2925SMatthew Knepley 
1707bee2925SMatthew Knepley         id   = EG_indexBodyTopo(body, edge);CHKERRQ(ierr);
171266cfabeSMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "            EDGE ID: %d (%d)\n", id, e);CHKERRQ(ierr);
1727bee2925SMatthew Knepley 
173266cfabeSMatthew G. Knepley         ierr = EG_getRange(edge, range, &peri);CHKERRQ(ierr);
1747bee2925SMatthew Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "  Range = %lf, %lf, %lf, %lf \n", range[0], range[1], range[2], range[3]);
1757bee2925SMatthew Knepley 
1767bee2925SMatthew Knepley         /* Get NODE info which associated with the current EDGE */
1777bee2925SMatthew Knepley         ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);CHKERRQ(ierr);
178266cfabeSMatthew G. Knepley         if (mtype == DEGENERATE) {
179266cfabeSMatthew G. Knepley           ierr = PetscPrintf(PETSC_COMM_SELF, "  EDGE %d is DEGENERATE \n", id);CHKERRQ(ierr);
180266cfabeSMatthew G. Knepley         } else {
181266cfabeSMatthew G. Knepley           params[0] = range[0];
182266cfabeSMatthew G. Knepley           ierr = EG_evaluate(edge, params, result);CHKERRQ(ierr);
183266cfabeSMatthew G. Knepley           ierr = PetscPrintf(PETSC_COMM_SELF, "   between (%lf, %lf, %lf)", result[0], result[1], result[2]);
184266cfabeSMatthew G. Knepley           params[0] = range[1];
185266cfabeSMatthew G. Knepley           ierr = EG_evaluate(edge, params, result);CHKERRQ(ierr);
186266cfabeSMatthew G. Knepley           ierr = PetscPrintf(PETSC_COMM_SELF, " and (%lf, %lf, %lf)\n", result[0], result[1], result[2]);
187266cfabeSMatthew G. Knepley         }
1887bee2925SMatthew Knepley 
1897bee2925SMatthew Knepley         for (v = 0; v < Nv; ++v) {
1907bee2925SMatthew Knepley           ego    vertex = nobjs[v];
1917bee2925SMatthew Knepley           double limits[4];
1927bee2925SMatthew Knepley           int    dummy;
1937bee2925SMatthew Knepley 
1947bee2925SMatthew Knepley           ierr = EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses);CHKERRQ(ierr);
1957bee2925SMatthew Knepley           id   = EG_indexBodyTopo(body, vertex);
1967bee2925SMatthew Knepley           ierr = PetscPrintf(PETSC_COMM_SELF, "              NODE ID: %d \n", id);CHKERRQ(ierr);
1977bee2925SMatthew Knepley           ierr = PetscPrintf(PETSC_COMM_SELF, "                 (x, y, z) = (%lf, %lf, %lf) \n", limits[0], limits[1], limits[2]);
1987bee2925SMatthew Knepley 
1997bee2925SMatthew Knepley           point[0] = point[0] + limits[0];
2007bee2925SMatthew Knepley           point[1] = point[1] + limits[1];
2017bee2925SMatthew Knepley           point[2] = point[2] + limits[2];
2027bee2925SMatthew Knepley         }
2037bee2925SMatthew Knepley       }
2047bee2925SMatthew Knepley     }
2057bee2925SMatthew Knepley     EG_free(lobjs);
2067bee2925SMatthew Knepley   }
2077bee2925SMatthew Knepley   PetscFunctionReturn(0);
2087bee2925SMatthew Knepley }
2097bee2925SMatthew Knepley 
2107bee2925SMatthew Knepley static PetscErrorCode DMPlexEGADSDestroy_Private(void *context)
2117bee2925SMatthew Knepley {
2127bee2925SMatthew Knepley   if (context) EG_close((ego) context);
2137bee2925SMatthew Knepley   return 0;
2147bee2925SMatthew Knepley }
2157bee2925SMatthew Knepley 
2167bee2925SMatthew Knepley static PetscErrorCode DMPlexCreateEGADS(MPI_Comm comm, ego context, ego model, DM *newdm)
2177bee2925SMatthew Knepley {
2187bee2925SMatthew Knepley   DMLabel        bodyLabel, faceLabel, edgeLabel;
2197bee2925SMatthew Knepley   PetscInt       cStart, cEnd, c;
2207bee2925SMatthew Knepley   /* EGADSLite variables */
2217bee2925SMatthew Knepley   ego            geom, *bodies, *objs, *nobjs, *mobjs, *lobjs;
2227bee2925SMatthew Knepley   int            oclass, mtype, nbodies, *senses;
2237bee2925SMatthew Knepley   int            b;
2247bee2925SMatthew Knepley   /* PETSc variables */
2257bee2925SMatthew Knepley   DM             dm;
2267bee2925SMatthew Knepley   PetscHMapI     edgeMap = NULL;
2277bee2925SMatthew Knepley   PetscInt       dim = -1, cdim = -1, numCorners = 0, maxCorners = 0, numVertices = 0, newVertices = 0, numEdges = 0, numCells = 0, newCells = 0, numQuads = 0, cOff = 0, fOff = 0;
2287bee2925SMatthew Knepley   PetscInt      *cells  = NULL, *cone = NULL;
2297bee2925SMatthew Knepley   PetscReal     *coords = NULL;
2307bee2925SMatthew Knepley   PetscMPIInt    rank;
2317bee2925SMatthew Knepley   PetscErrorCode ierr;
2327bee2925SMatthew Knepley 
2337bee2925SMatthew Knepley   PetscFunctionBegin;
2347bee2925SMatthew Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2357bee2925SMatthew Knepley   if (!rank) {
236266cfabeSMatthew G. Knepley     const PetscInt debug = 0;
237266cfabeSMatthew G. Knepley 
2387bee2925SMatthew Knepley     /* ---------------------------------------------------------------------------------------------------
2397bee2925SMatthew Knepley     Generate Petsc Plex
2407bee2925SMatthew Knepley       Get all Nodes in model, record coordinates in a correctly formatted array
2417bee2925SMatthew Knepley       Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array
2427bee2925SMatthew Knepley       We need to uniformly refine the initial geometry to guarantee a valid mesh
2437bee2925SMatthew Knepley     */
2447bee2925SMatthew Knepley 
2457bee2925SMatthew Knepley     /* Calculate cell and vertex sizes */
2467bee2925SMatthew Knepley     ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses);CHKERRQ(ierr);
2477bee2925SMatthew Knepley     ierr = PetscHMapICreate(&edgeMap);CHKERRQ(ierr);
2487bee2925SMatthew Knepley     numEdges = 0;
2497bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
2507bee2925SMatthew Knepley       ego body = bodies[b];
2517bee2925SMatthew Knepley       int id, Nl, l, Nv, v;
2527bee2925SMatthew Knepley 
2537bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
2547bee2925SMatthew Knepley       for (l = 0; l < Nl; ++l) {
2557bee2925SMatthew Knepley         ego loop = lobjs[l];
2567bee2925SMatthew Knepley         int Ner  = 0, Ne, e, Nc;
2577bee2925SMatthew Knepley 
2587bee2925SMatthew Knepley         ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
2597bee2925SMatthew Knepley         for (e = 0; e < Ne; ++e) {
2607bee2925SMatthew Knepley           ego edge = objs[e];
2617bee2925SMatthew Knepley           int Nv, id;
2627bee2925SMatthew Knepley           PetscHashIter iter;
2637bee2925SMatthew Knepley           PetscBool     found;
2647bee2925SMatthew Knepley 
2657bee2925SMatthew Knepley           ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);
2667bee2925SMatthew Knepley           if (mtype == DEGENERATE) continue;
2677bee2925SMatthew Knepley           id   = EG_indexBodyTopo(body, edge);CHKERRQ(ierr);
2687bee2925SMatthew Knepley           ierr = PetscHMapIFind(edgeMap, id-1, &iter, &found);CHKERRQ(ierr);
2697bee2925SMatthew Knepley           if (!found) {ierr = PetscHMapISet(edgeMap, id-1, numEdges++);CHKERRQ(ierr);}
2707bee2925SMatthew Knepley           ++Ner;
2717bee2925SMatthew Knepley         }
2727bee2925SMatthew Knepley         if (Ner == 2)      {Nc = 2;}
2737bee2925SMatthew Knepley         else if (Ner == 3) {Nc = 4;}
2747bee2925SMatthew Knepley         else if (Ner == 4) {Nc = 8; ++numQuads;}
2757bee2925SMatthew Knepley         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot support loop with %d edges", Ner);
2767bee2925SMatthew Knepley         numCells += Nc;
2777bee2925SMatthew Knepley         newCells += Nc-1;
2787bee2925SMatthew Knepley         maxCorners = PetscMax(Ner*2+1, maxCorners);
2797bee2925SMatthew Knepley       }
2807bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs);CHKERRQ(ierr);
2817bee2925SMatthew Knepley       for (v = 0; v < Nv; ++v) {
2827bee2925SMatthew Knepley         ego vertex = nobjs[v];
2837bee2925SMatthew Knepley 
2847bee2925SMatthew Knepley         id = EG_indexBodyTopo(body, vertex);
2857bee2925SMatthew Knepley         /* TODO: Instead of assuming contiguous ids, we could use a hash table */
2867bee2925SMatthew Knepley         numVertices = PetscMax(id, numVertices);
2877bee2925SMatthew Knepley       }
2887bee2925SMatthew Knepley       EG_free(lobjs);
2897bee2925SMatthew Knepley       EG_free(nobjs);
2907bee2925SMatthew Knepley     }
2917bee2925SMatthew Knepley     ierr = PetscHMapIGetSize(edgeMap, &numEdges);CHKERRQ(ierr);
2927bee2925SMatthew Knepley     newVertices  = numEdges + numQuads;
2937bee2925SMatthew Knepley     numVertices += newVertices;
2947bee2925SMatthew Knepley 
2957bee2925SMatthew Knepley     dim        = 2; /* Assume 3D Models :: Need to update to handle 2D Models in the future */
2967bee2925SMatthew Knepley     cdim       = 3; /* Assume 3D Models :: Need to update to handle 2D Models in the future */
2977bee2925SMatthew Knepley     numCorners = 3; /* Split cells into triangles */
2987bee2925SMatthew Knepley     ierr = PetscMalloc3(numVertices*cdim, &coords, numCells*numCorners, &cells, maxCorners, &cone);CHKERRQ(ierr);
2997bee2925SMatthew Knepley 
3007bee2925SMatthew Knepley     /* Get vertex coordinates */
3017bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
3027bee2925SMatthew Knepley       ego body = bodies[b];
3037bee2925SMatthew Knepley       int id, Nv, v;
3047bee2925SMatthew Knepley 
3057bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs);CHKERRQ(ierr);
3067bee2925SMatthew Knepley       for (v = 0; v < Nv; ++v) {
3077bee2925SMatthew Knepley         ego    vertex = nobjs[v];
3087bee2925SMatthew Knepley         double limits[4];
3097bee2925SMatthew Knepley         int    dummy;
3107bee2925SMatthew Knepley 
3117bee2925SMatthew Knepley         ierr = EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses);CHKERRQ(ierr);
3127bee2925SMatthew Knepley         id   = EG_indexBodyTopo(body, vertex);CHKERRQ(ierr);
3137bee2925SMatthew Knepley         coords[(id-1)*cdim+0] = limits[0];
3147bee2925SMatthew Knepley         coords[(id-1)*cdim+1] = limits[1];
3157bee2925SMatthew Knepley         coords[(id-1)*cdim+2] = limits[2];
3167bee2925SMatthew Knepley       }
3177bee2925SMatthew Knepley       EG_free(nobjs);
3187bee2925SMatthew Knepley     }
3197bee2925SMatthew Knepley     ierr = PetscHMapIClear(edgeMap);CHKERRQ(ierr);
3207bee2925SMatthew Knepley     fOff     = numVertices - newVertices + numEdges;
3217bee2925SMatthew Knepley     numEdges = 0;
3227bee2925SMatthew Knepley     numQuads = 0;
3237bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
3247bee2925SMatthew Knepley       ego body = bodies[b];
3257bee2925SMatthew Knepley       int Nl, l;
3267bee2925SMatthew Knepley 
3277bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
3287bee2925SMatthew Knepley       for (l = 0; l < Nl; ++l) {
3297bee2925SMatthew Knepley         ego loop = lobjs[l];
3307bee2925SMatthew Knepley         int lid, Ner = 0, Ne, e;
3317bee2925SMatthew Knepley 
3327bee2925SMatthew Knepley         lid  = EG_indexBodyTopo(body, loop);CHKERRQ(ierr);
3337bee2925SMatthew Knepley         ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
3347bee2925SMatthew Knepley         for (e = 0; e < Ne; ++e) {
3357bee2925SMatthew Knepley           ego       edge = objs[e];
3367bee2925SMatthew Knepley           int       eid, Nv;
3377bee2925SMatthew Knepley           PetscHashIter iter;
3387bee2925SMatthew Knepley           PetscBool     found;
3397bee2925SMatthew Knepley 
3407bee2925SMatthew Knepley           ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);
3417bee2925SMatthew Knepley           if (mtype == DEGENERATE) continue;
3427bee2925SMatthew Knepley           ++Ner;
3437bee2925SMatthew Knepley           eid  = EG_indexBodyTopo(body, edge);CHKERRQ(ierr);
3447bee2925SMatthew Knepley           ierr = PetscHMapIFind(edgeMap, eid-1, &iter, &found);CHKERRQ(ierr);
3457bee2925SMatthew Knepley           if (!found) {
3467bee2925SMatthew Knepley             PetscInt v = numVertices - newVertices + numEdges;
3477bee2925SMatthew Knepley             double range[4], params[3] = {0., 0., 0.}, result[18];
3487bee2925SMatthew Knepley             int    periodic[2];
3497bee2925SMatthew Knepley 
3507bee2925SMatthew Knepley             ierr = PetscHMapISet(edgeMap, eid-1, numEdges++);CHKERRQ(ierr);
3517bee2925SMatthew Knepley             ierr = EG_getRange(edge, range, periodic);CHKERRQ(ierr);
3527bee2925SMatthew Knepley             params[0] = 0.5*(range[0] + range[1]);
3537bee2925SMatthew Knepley             ierr = EG_evaluate(edge, params, result);CHKERRQ(ierr);
3547bee2925SMatthew Knepley             coords[v*cdim+0] = result[0];
3557bee2925SMatthew Knepley             coords[v*cdim+1] = result[1];
3567bee2925SMatthew Knepley             coords[v*cdim+2] = result[2];
3577bee2925SMatthew Knepley           }
3587bee2925SMatthew Knepley         }
3597bee2925SMatthew Knepley         if (Ner == 4) {
3607bee2925SMatthew Knepley           PetscInt v = fOff + numQuads++;
361266cfabeSMatthew G. Knepley           ego     *fobjs, face;
3627bee2925SMatthew Knepley           double   range[4], params[3] = {0., 0., 0.}, result[18];
363266cfabeSMatthew G. Knepley           int      Nf, fid, periodic[2];
3647bee2925SMatthew Knepley 
3657bee2925SMatthew Knepley           ierr = EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs);CHKERRQ(ierr);
366266cfabeSMatthew G. Knepley           face = fobjs[0];
367266cfabeSMatthew G. Knepley           fid  = EG_indexBodyTopo(body, face);CHKERRQ(ierr);
368266cfabeSMatthew G. Knepley           if (Nf != 1) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Loop %d has %d faces, instead of 1 (%d)", lid-1, Nf, fid);
369266cfabeSMatthew G. Knepley           ierr = EG_getRange(face, range, periodic);CHKERRQ(ierr);
3707bee2925SMatthew Knepley           params[0] = 0.5*(range[0] + range[1]);
3717bee2925SMatthew Knepley           params[1] = 0.5*(range[2] + range[3]);
372266cfabeSMatthew G. Knepley           ierr = EG_evaluate(face, params, result);CHKERRQ(ierr);
3737bee2925SMatthew Knepley           coords[v*cdim+0] = result[0];
3747bee2925SMatthew Knepley           coords[v*cdim+1] = result[1];
3757bee2925SMatthew Knepley           coords[v*cdim+2] = result[2];
3767bee2925SMatthew Knepley         }
3777bee2925SMatthew Knepley       }
3787bee2925SMatthew Knepley     }
3797bee2925SMatthew Knepley     if (numEdges + numQuads != newVertices) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of new vertices %D != %D previous count", numEdges + numQuads, newVertices);
3807bee2925SMatthew Knepley     CHKMEMQ;
3817bee2925SMatthew Knepley 
3827bee2925SMatthew Knepley     /* Get cell vertices by traversing loops */
3837bee2925SMatthew Knepley     numQuads = 0;
3847bee2925SMatthew Knepley     cOff     = 0;
3857bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
3867bee2925SMatthew Knepley       ego body = bodies[b];
3877bee2925SMatthew Knepley       int id, Nl, l;
3887bee2925SMatthew Knepley 
3897bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
3907bee2925SMatthew Knepley       for (l = 0; l < Nl; ++l) {
3917bee2925SMatthew Knepley         ego loop = lobjs[l];
3927bee2925SMatthew Knepley         int lid, Ner = 0, Ne, e, nc = 0, c, Nt, t;
3937bee2925SMatthew Knepley 
3947bee2925SMatthew Knepley         lid  = EG_indexBodyTopo(body, loop);CHKERRQ(ierr);
3957bee2925SMatthew Knepley         ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
3967bee2925SMatthew Knepley 
3977bee2925SMatthew Knepley         for (e = 0; e < Ne; ++e) {
3987bee2925SMatthew Knepley           ego edge = objs[e];
3997bee2925SMatthew Knepley           int points[3];
4007bee2925SMatthew Knepley           int eid, Nv, v, tmp;
4017bee2925SMatthew Knepley 
4027bee2925SMatthew Knepley           eid  = EG_indexBodyTopo(body, edge);
4037bee2925SMatthew Knepley           ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);
404266cfabeSMatthew G. Knepley           if (mtype == DEGENERATE) continue;
405266cfabeSMatthew G. Knepley           else                     ++Ner;
4067bee2925SMatthew Knepley           if (Nv != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Edge %d has %d vertices != 2", eid, Nv);
4077bee2925SMatthew Knepley 
4087bee2925SMatthew Knepley           for (v = 0; v < Nv; ++v) {
4097bee2925SMatthew Knepley             ego vertex = nobjs[v];
4107bee2925SMatthew Knepley 
4117bee2925SMatthew Knepley             id = EG_indexBodyTopo(body, vertex);
4127bee2925SMatthew Knepley             points[v*2] = id-1;
4137bee2925SMatthew Knepley           }
4147bee2925SMatthew Knepley           {
4157bee2925SMatthew Knepley             PetscInt edgeNum;
4167bee2925SMatthew Knepley 
4177bee2925SMatthew Knepley             ierr = PetscHMapIGet(edgeMap, eid-1, &edgeNum);CHKERRQ(ierr);
4187bee2925SMatthew Knepley             points[1] = numVertices - newVertices + edgeNum;
4197bee2925SMatthew Knepley           }
4207bee2925SMatthew Knepley           /* EGADS loops are not oriented, but seem to be in order, so we must piece them together */
4217bee2925SMatthew Knepley           if (!nc) {
4227bee2925SMatthew Knepley             for (v = 0; v < Nv+1; ++v) cone[nc++] = points[v];
4237bee2925SMatthew Knepley           } else {
4247bee2925SMatthew Knepley             if (cone[nc-1] == points[0])      {cone[nc++] = points[1]; if (cone[0] != points[2]) cone[nc++] = points[2];}
4257bee2925SMatthew Knepley             else if (cone[nc-1] == points[2]) {cone[nc++] = points[1]; if (cone[0] != points[0]) cone[nc++] = points[0];}
4267bee2925SMatthew Knepley             else if (cone[nc-3] == points[0]) {tmp = cone[nc-3]; cone[nc-3] = cone[nc-1]; cone[nc-1] = tmp; cone[nc++] = points[1]; if (cone[0] != points[2]) cone[nc++] = points[2];}
4277bee2925SMatthew Knepley             else if (cone[nc-3] == points[2]) {tmp = cone[nc-3]; cone[nc-3] = cone[nc-1]; cone[nc-1] = tmp; cone[nc++] = points[1]; if (cone[0] != points[0]) cone[nc++] = points[0];}
4287bee2925SMatthew Knepley             else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %d does not match its predecessor", eid);
4297bee2925SMatthew Knepley           }
4307bee2925SMatthew Knepley         }
4317bee2925SMatthew Knepley         if (nc != 2*Ner)     SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %D != %D", nc, 2*Ner);
4327bee2925SMatthew Knepley         if (Ner == 4) {cone[nc++] = numVertices - newVertices + numEdges + numQuads++;}
4337bee2925SMatthew Knepley         if (nc > maxCorners) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %D > %D max", nc, maxCorners);
4347bee2925SMatthew Knepley         /* Triangulate the loop */
4357bee2925SMatthew Knepley         switch (Ner) {
4367bee2925SMatthew Knepley           case 2: /* Bi-Segment -> 2 triangles */
4377bee2925SMatthew Knepley             Nt = 2;
4387bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[0];
4397bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[1];
4407bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[2];
4417bee2925SMatthew Knepley             ++cOff;
4427bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[0];
4437bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[2];
4447bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[3];
4457bee2925SMatthew Knepley             ++cOff;
4467bee2925SMatthew Knepley             break;
4477bee2925SMatthew Knepley           case 3: /* Triangle   -> 4 triangles */
4487bee2925SMatthew Knepley             Nt = 4;
4497bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[0];
4507bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[1];
4517bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[5];
4527bee2925SMatthew Knepley             ++cOff;
4537bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[1];
4547bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[2];
4557bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[3];
4567bee2925SMatthew Knepley             ++cOff;
4577bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[5];
4587bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[3];
4597bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[4];
4607bee2925SMatthew Knepley             ++cOff;
4617bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[1];
4627bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[3];
4637bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[5];
4647bee2925SMatthew Knepley             ++cOff;
4657bee2925SMatthew Knepley             break;
4667bee2925SMatthew Knepley           case 4: /* Quad       -> 8 triangles */
4677bee2925SMatthew Knepley             Nt = 8;
4687bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[0];
4697bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[1];
4707bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[7];
4717bee2925SMatthew Knepley             ++cOff;
4727bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[1];
4737bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[2];
4747bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[3];
4757bee2925SMatthew Knepley             ++cOff;
4767bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[3];
4777bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[4];
4787bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[5];
4797bee2925SMatthew Knepley             ++cOff;
4807bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[5];
4817bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[6];
4827bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[7];
4837bee2925SMatthew Knepley             ++cOff;
4847bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[8];
4857bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[1];
4867bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[3];
4877bee2925SMatthew Knepley             ++cOff;
4887bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[8];
4897bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[3];
4907bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[5];
4917bee2925SMatthew Knepley             ++cOff;
4927bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[8];
4937bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[5];
4947bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[7];
4957bee2925SMatthew Knepley             ++cOff;
4967bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[8];
4977bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[7];
4987bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[1];
4997bee2925SMatthew Knepley             ++cOff;
5007bee2925SMatthew Knepley             break;
5017bee2925SMatthew Knepley           default: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d edges, which we do not support", lid, Ner);
5027bee2925SMatthew Knepley         }
503266cfabeSMatthew G. Knepley         if (debug) {
5047bee2925SMatthew Knepley           for (t = 0; t < Nt; ++t) {
5057bee2925SMatthew Knepley             ierr = PetscPrintf(PETSC_COMM_SELF, "  LOOP Corner NODEs Triangle %D (", t);CHKERRQ(ierr);
5067bee2925SMatthew Knepley             for (c = 0; c < numCorners; ++c) {
5077bee2925SMatthew Knepley               if (c > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, ", ");CHKERRQ(ierr);}
5087bee2925SMatthew Knepley               ierr = PetscPrintf(PETSC_COMM_SELF, "%D", cells[(cOff-Nt+t)*numCorners+c]);CHKERRQ(ierr);
5097bee2925SMatthew Knepley             }
5107bee2925SMatthew Knepley             ierr = PetscPrintf(PETSC_COMM_SELF, ")\n");CHKERRQ(ierr);
5117bee2925SMatthew Knepley           }
5127bee2925SMatthew Knepley         }
513266cfabeSMatthew G. Knepley       }
5147bee2925SMatthew Knepley       EG_free(lobjs);
5157bee2925SMatthew Knepley     }
5167bee2925SMatthew Knepley   }
5177bee2925SMatthew Knepley   if (cOff != numCells) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count of total cells %D != %D previous count", cOff, numCells);
5187bee2925SMatthew Knepley   ierr = DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numVertices, numCorners, PETSC_TRUE, cells, cdim, coords, &dm);CHKERRQ(ierr);
5197bee2925SMatthew Knepley   ierr = PetscFree3(coords, cells, cone);CHKERRQ(ierr);
520266cfabeSMatthew G. Knepley   ierr = PetscInfo2(dm, " Total Number of Unique Cells    = %D (%D)\n", numCells, newCells);CHKERRQ(ierr);
521266cfabeSMatthew G. Knepley   ierr = PetscInfo2(dm, " Total Number of Unique Vertices = %D (%D)\n", numVertices, newVertices);CHKERRQ(ierr);
5227bee2925SMatthew Knepley   /* Embed EGADS model in DM */
5237bee2925SMatthew Knepley   {
5247bee2925SMatthew Knepley     PetscContainer modelObj, contextObj;
5257bee2925SMatthew Knepley 
5267bee2925SMatthew Knepley     ierr = PetscContainerCreate(PETSC_COMM_SELF, &modelObj);CHKERRQ(ierr);
5277bee2925SMatthew Knepley     ierr = PetscContainerSetPointer(modelObj, model);CHKERRQ(ierr);
5287bee2925SMatthew Knepley     ierr = PetscObjectCompose((PetscObject) dm, "EGADS Model", (PetscObject) modelObj);CHKERRQ(ierr);
5297bee2925SMatthew Knepley     ierr = PetscContainerDestroy(&modelObj);CHKERRQ(ierr);
5307bee2925SMatthew Knepley 
5317bee2925SMatthew Knepley     ierr = PetscContainerCreate(PETSC_COMM_SELF, &contextObj);CHKERRQ(ierr);
5327bee2925SMatthew Knepley     ierr = PetscContainerSetPointer(contextObj, context);CHKERRQ(ierr);
5337bee2925SMatthew Knepley     ierr = PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private);CHKERRQ(ierr);
5347bee2925SMatthew Knepley     ierr = PetscObjectCompose((PetscObject) dm, "EGADS Context", (PetscObject) contextObj);CHKERRQ(ierr);
5357bee2925SMatthew Knepley     ierr = PetscContainerDestroy(&contextObj);CHKERRQ(ierr);
5367bee2925SMatthew Knepley   }
5377bee2925SMatthew Knepley   /* Label points */
5387bee2925SMatthew Knepley   ierr = DMCreateLabel(dm, "EGADS Body ID");CHKERRQ(ierr);
5397bee2925SMatthew Knepley   ierr = DMGetLabel(dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
5407bee2925SMatthew Knepley   ierr = DMCreateLabel(dm, "EGADS Face ID");CHKERRQ(ierr);
5417bee2925SMatthew Knepley   ierr = DMGetLabel(dm, "EGADS Face ID", &faceLabel);CHKERRQ(ierr);
5427bee2925SMatthew Knepley   ierr = DMCreateLabel(dm, "EGADS Edge ID");CHKERRQ(ierr);
5437bee2925SMatthew Knepley   ierr = DMGetLabel(dm, "EGADS Edge ID", &edgeLabel);CHKERRQ(ierr);
5447bee2925SMatthew Knepley   cOff = 0;
5457bee2925SMatthew Knepley   for (b = 0; b < nbodies; ++b) {
5467bee2925SMatthew Knepley     ego body = bodies[b];
5477bee2925SMatthew Knepley     int id, Nl, l;
5487bee2925SMatthew Knepley 
5497bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
5507bee2925SMatthew Knepley     for (l = 0; l < Nl; ++l) {
5517bee2925SMatthew Knepley       ego  loop = lobjs[l];
5527bee2925SMatthew Knepley       ego *fobjs;
5537bee2925SMatthew Knepley       int  lid, Nf, fid, Ner = 0, Ne, e, Nt = 0, t;
5547bee2925SMatthew Knepley 
5557bee2925SMatthew Knepley       lid  = EG_indexBodyTopo(body, loop);CHKERRQ(ierr);
5567bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs);CHKERRQ(ierr);
5577bee2925SMatthew Knepley       if (Nf > 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf);CHKERRQ(ierr);
5587bee2925SMatthew Knepley       fid  = EG_indexBodyTopo(body, fobjs[0]);CHKERRQ(ierr);
5597bee2925SMatthew Knepley       EG_free(fobjs);
5607bee2925SMatthew Knepley       ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
5617bee2925SMatthew Knepley       for (e = 0; e < Ne; ++e) {
5627bee2925SMatthew Knepley         ego             edge = objs[e];
5637bee2925SMatthew Knepley         int             eid, Nv, v;
5647bee2925SMatthew Knepley         PetscInt        points[3], support[2], numEdges, edgeNum;
5657bee2925SMatthew Knepley         const PetscInt *edges;
5667bee2925SMatthew Knepley 
5677bee2925SMatthew Knepley         eid  = EG_indexBodyTopo(body, edge);
5687bee2925SMatthew Knepley         ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);
5697bee2925SMatthew Knepley         if (mtype == DEGENERATE) continue;
5707bee2925SMatthew Knepley         else                     ++Ner;
5717bee2925SMatthew Knepley         for (v = 0; v < Nv; ++v) {
5727bee2925SMatthew Knepley           ego vertex = nobjs[v];
5737bee2925SMatthew Knepley 
5747bee2925SMatthew Knepley           id   = EG_indexBodyTopo(body, vertex);
5757bee2925SMatthew Knepley           ierr = DMLabelSetValue(edgeLabel, numCells + id-1, eid);CHKERRQ(ierr);
5767bee2925SMatthew Knepley           points[v*2] = numCells + id-1;
5777bee2925SMatthew Knepley         }
5787bee2925SMatthew Knepley         ierr = PetscHMapIGet(edgeMap, eid-1, &edgeNum);CHKERRQ(ierr);
5797bee2925SMatthew Knepley         points[1] = numCells + numVertices - newVertices + edgeNum;
5807bee2925SMatthew Knepley 
5817bee2925SMatthew Knepley         ierr = DMLabelSetValue(edgeLabel, points[1], eid);CHKERRQ(ierr);
5827bee2925SMatthew Knepley         support[0] = points[0];
5837bee2925SMatthew Knepley         support[1] = points[1];
5847bee2925SMatthew Knepley         ierr = DMPlexGetJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr);
5857bee2925SMatthew Knepley         if (numEdges != 1) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vertices (%D, %D) should only bound 1 edge, not %D", support[0], support[1], numEdges);
5867bee2925SMatthew Knepley         ierr = DMLabelSetValue(edgeLabel, edges[0], eid);CHKERRQ(ierr);
5877bee2925SMatthew Knepley         ierr = DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr);
5887bee2925SMatthew Knepley         support[0] = points[1];
5897bee2925SMatthew Knepley         support[1] = points[2];
5907bee2925SMatthew Knepley         ierr = DMPlexGetJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr);
5917bee2925SMatthew Knepley         if (numEdges != 1) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vertices (%D, %D) should only bound 1 edge, not %D", support[0], support[1], numEdges);
5927bee2925SMatthew Knepley         ierr = DMLabelSetValue(edgeLabel, edges[0], eid);CHKERRQ(ierr);
5937bee2925SMatthew Knepley         ierr = DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr);
5947bee2925SMatthew Knepley       }
5957bee2925SMatthew Knepley       switch (Ner) {
5967bee2925SMatthew Knepley         case 2: Nt = 2;break;
5977bee2925SMatthew Knepley         case 3: Nt = 4;break;
5987bee2925SMatthew Knepley         case 4: Nt = 8;break;
5997bee2925SMatthew Knepley         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Loop with %d edges is unsupported", Ner);
6007bee2925SMatthew Knepley       }
6017bee2925SMatthew Knepley       for (t = 0; t < Nt; ++t) {
6027bee2925SMatthew Knepley         ierr = DMLabelSetValue(bodyLabel, cOff+t, b);CHKERRQ(ierr);
6037bee2925SMatthew Knepley         ierr = DMLabelSetValue(faceLabel, cOff+t, fid);CHKERRQ(ierr);
6047bee2925SMatthew Knepley       }
6057bee2925SMatthew Knepley       cOff += Nt;
6067bee2925SMatthew Knepley     }
6077bee2925SMatthew Knepley     EG_free(lobjs);
6087bee2925SMatthew Knepley   }
6097bee2925SMatthew Knepley   ierr = PetscHMapIDestroy(&edgeMap);CHKERRQ(ierr);
6107bee2925SMatthew Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
6117bee2925SMatthew Knepley   for (c = cStart; c < cEnd; ++c) {
6127bee2925SMatthew Knepley     PetscInt *closure = NULL;
6137bee2925SMatthew Knepley     PetscInt  clSize, cl, bval, fval;
6147bee2925SMatthew Knepley 
6157bee2925SMatthew Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
6167bee2925SMatthew Knepley     ierr = DMLabelGetValue(bodyLabel, c, &bval);CHKERRQ(ierr);
6177bee2925SMatthew Knepley     ierr = DMLabelGetValue(faceLabel, c, &fval);CHKERRQ(ierr);
6187bee2925SMatthew Knepley     for (cl = 0; cl < clSize*2; cl += 2) {
6197bee2925SMatthew Knepley       ierr = DMLabelSetValue(bodyLabel, closure[cl], bval);CHKERRQ(ierr);
6207bee2925SMatthew Knepley       ierr = DMLabelSetValue(faceLabel, closure[cl], fval);CHKERRQ(ierr);
6217bee2925SMatthew Knepley     }
6227bee2925SMatthew Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
6237bee2925SMatthew Knepley   }
6247bee2925SMatthew Knepley   *newdm = dm;
6257bee2925SMatthew Knepley   PetscFunctionReturn(0);
6267bee2925SMatthew Knepley }
6277bee2925SMatthew Knepley #endif
6287bee2925SMatthew Knepley 
6297bee2925SMatthew Knepley /*@C
6307bee2925SMatthew Knepley   DMPlexCreateEGADSFromFile - Create a DMPlex mesh from an EGADSLite file.
6317bee2925SMatthew Knepley 
6327bee2925SMatthew Knepley   Collective
6337bee2925SMatthew Knepley 
6347bee2925SMatthew Knepley   Input Parameters:
6357bee2925SMatthew Knepley + comm     - The MPI communicator
6367bee2925SMatthew Knepley - filename - The name of the ExodusII file
6377bee2925SMatthew Knepley 
6387bee2925SMatthew Knepley   Output Parameter:
6397bee2925SMatthew Knepley . dm       - The DM object representing the mesh
6407bee2925SMatthew Knepley 
6417bee2925SMatthew Knepley   Level: beginner
6427bee2925SMatthew Knepley 
6437bee2925SMatthew Knepley .seealso: DMPLEX, DMCreate(), DMPlexCreateEGADS()
6447bee2925SMatthew Knepley @*/
6457bee2925SMatthew Knepley PetscErrorCode DMPlexCreateEGADSFromFile(MPI_Comm comm, const char filename[], DM *dm)
6467bee2925SMatthew Knepley {
6477bee2925SMatthew Knepley   PetscMPIInt    rank;
6487bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS)
6497bee2925SMatthew Knepley   ego            context= NULL, model = NULL;
6507bee2925SMatthew Knepley #endif
651266cfabeSMatthew G. Knepley   PetscBool      printModel = PETSC_FALSE;
6527bee2925SMatthew Knepley   PetscErrorCode ierr;
6537bee2925SMatthew Knepley 
6547bee2925SMatthew Knepley   PetscFunctionBegin;
6557bee2925SMatthew Knepley   PetscValidCharPointer(filename, 2);
6567bee2925SMatthew Knepley   ierr = PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_print_model", &printModel, NULL);CHKERRQ(ierr);
6577bee2925SMatthew Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
6587bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS)
6597bee2925SMatthew Knepley   if (!rank) {
6607bee2925SMatthew Knepley 
6617bee2925SMatthew Knepley     ierr = EG_open(&context);CHKERRQ(ierr);
6627bee2925SMatthew Knepley     ierr = EG_loadModel(context, 0, filename, &model);CHKERRQ(ierr);
6637bee2925SMatthew Knepley     if (printModel) {ierr = DMPlexEGADSPrintModel(model);CHKERRQ(ierr);}
6647bee2925SMatthew Knepley 
6657bee2925SMatthew Knepley   }
6667bee2925SMatthew Knepley   ierr = DMPlexCreateEGADS(comm, context, model, dm);CHKERRQ(ierr);
667*c03d1177SJose E. Roman   PetscFunctionReturn(0);
6687bee2925SMatthew Knepley #else
6697bee2925SMatthew Knepley   SETERRQ(comm, PETSC_ERR_SUP, "This method requires EGADSLite support. Reconfigure using --download-egads");
6707bee2925SMatthew Knepley #endif
6717bee2925SMatthew Knepley }
672