xref: /petsc/src/dm/impls/plex/plexegads.c (revision d410b0cf18e1798d3d4c14858e0c2ffdbe2fea69)
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 
21c1cad2e7SMatthew G. Knepley #ifdef PETSC_HAVE_EGADS
22c1cad2e7SMatthew G. Knepley PETSC_INTERN PetscErrorCode DMPlexSnapToGeomModel_EGADS_Internal(DM, PetscInt, ego, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
23c1cad2e7SMatthew G. Knepley PETSC_INTERN PetscErrorCode DMPlexSnapToGeomModel_EGADSLite_Internal(DM, PetscInt, ego, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]);
24c1cad2e7SMatthew G. Knepley 
25c1cad2e7SMatthew G. Knepley PetscErrorCode DMPlexSnapToGeomModel_EGADS_Internal(DM dm, PetscInt p, ego model, PetscInt bodyID, PetscInt faceID, PetscInt edgeID, const PetscScalar mcoords[], PetscScalar gcoords[])
26c1cad2e7SMatthew G. Knepley {
27c1cad2e7SMatthew G. Knepley   DM             cdm;
28c1cad2e7SMatthew G. Knepley   ego           *bodies;
29c1cad2e7SMatthew G. Knepley   ego            geom, body, obj;
30a5b23f4aSJose E. Roman   /* result has to hold derivatives, along with the value */
31c1cad2e7SMatthew G. Knepley   double         params[3], result[18], paramsV[16*3], resultV[16*3], range[4];
32c1cad2e7SMatthew G. Knepley   int            Nb, oclass, mtype, *senses, peri;
33c1cad2e7SMatthew G. Knepley   Vec            coordinatesLocal;
34c1cad2e7SMatthew G. Knepley   PetscScalar   *coords = NULL;
35c1cad2e7SMatthew G. Knepley   PetscInt       Nv, v, Np = 0, pm;
36c1cad2e7SMatthew G. Knepley   PetscInt       dE, d;
37c1cad2e7SMatthew G. Knepley   PetscErrorCode ierr;
38c1cad2e7SMatthew G. Knepley 
39c1cad2e7SMatthew G. Knepley   PetscFunctionBeginHot;
40c1cad2e7SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
41c1cad2e7SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &dE);CHKERRQ(ierr);
42c1cad2e7SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinatesLocal);CHKERRQ(ierr);
43c1cad2e7SMatthew G. Knepley   ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
44c1cad2e7SMatthew G. Knepley   if (bodyID >= Nb) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %D is not in [0, %d)", bodyID, Nb);
45c1cad2e7SMatthew G. Knepley   body = bodies[bodyID];
46c1cad2e7SMatthew G. Knepley 
47c1cad2e7SMatthew G. Knepley   if (edgeID >= 0)      {ierr = EG_objectBodyTopo(body, EDGE, edgeID, &obj);CHKERRQ(ierr); Np = 1;}
48c1cad2e7SMatthew G. Knepley   else if (faceID >= 0) {ierr = EG_objectBodyTopo(body, FACE, faceID, &obj);CHKERRQ(ierr); Np = 2;}
49c1cad2e7SMatthew G. Knepley   else {
50c1cad2e7SMatthew G. Knepley     for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
51c1cad2e7SMatthew G. Knepley     PetscFunctionReturn(0);
52c1cad2e7SMatthew G. Knepley   }
53c1cad2e7SMatthew G. Knepley 
54c1cad2e7SMatthew G. Knepley   /* Calculate parameters (t or u,v) for vertices */
55c1cad2e7SMatthew G. Knepley   ierr = DMPlexVecGetClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords);CHKERRQ(ierr);
56c1cad2e7SMatthew G. Knepley   Nv  /= dE;
57c1cad2e7SMatthew G. Knepley   if (Nv == 1) {
58c1cad2e7SMatthew G. Knepley     ierr = DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords);CHKERRQ(ierr);
59c1cad2e7SMatthew G. Knepley     for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
60c1cad2e7SMatthew G. Knepley     PetscFunctionReturn(0);
61c1cad2e7SMatthew G. Knepley   }
62c1cad2e7SMatthew G. Knepley   if (Nv > 16) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %D coordinates associated to point %D", Nv, p);
63c1cad2e7SMatthew G. Knepley 
64c1cad2e7SMatthew G. Knepley   /* Correct EGADSlite 2pi bug when calculating nearest point on Periodic Surfaces */
65c1cad2e7SMatthew G. Knepley   ierr = EG_getRange(obj, range, &peri);CHKERRQ(ierr);
66c1cad2e7SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
67c1cad2e7SMatthew G. Knepley     ierr = EG_invEvaluate(obj, &coords[v*dE], &paramsV[v*3], &resultV[v*3]);CHKERRQ(ierr);
68c1cad2e7SMatthew G. Knepley #if 1
69c1cad2e7SMatthew G. Knepley     if (peri > 0) {
70c1cad2e7SMatthew G. Knepley       if      (paramsV[v*3+0] + 1.e-4 < range[0]) {paramsV[v*3+0] += 2. * PETSC_PI;}
71c1cad2e7SMatthew G. Knepley       else if (paramsV[v*3+0] - 1.e-4 > range[1]) {paramsV[v*3+0] -= 2. * PETSC_PI;}
72c1cad2e7SMatthew G. Knepley     }
73c1cad2e7SMatthew G. Knepley     if (peri > 1) {
74c1cad2e7SMatthew G. Knepley       if      (paramsV[v*3+1] + 1.e-4 < range[2]) {paramsV[v*3+1] += 2. * PETSC_PI;}
75c1cad2e7SMatthew G. Knepley       else if (paramsV[v*3+1] - 1.e-4 > range[3]) {paramsV[v*3+1] -= 2. * PETSC_PI;}
76c1cad2e7SMatthew G. Knepley     }
77c1cad2e7SMatthew G. Knepley #endif
78c1cad2e7SMatthew G. Knepley   }
79c1cad2e7SMatthew G. Knepley   ierr = DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords);CHKERRQ(ierr);
80c1cad2e7SMatthew G. Knepley   /* Calculate parameters (t or u,v) for new vertex at edge midpoint */
81c1cad2e7SMatthew G. Knepley   for (pm = 0; pm < Np; ++pm) {
82c1cad2e7SMatthew G. Knepley     params[pm] = 0.;
83c1cad2e7SMatthew G. Knepley     for (v = 0; v < Nv; ++v) params[pm] += paramsV[v*3+pm];
84c1cad2e7SMatthew G. Knepley     params[pm] /= Nv;
85c1cad2e7SMatthew G. Knepley   }
86c1cad2e7SMatthew G. Knepley   if ((params[0] < range[0]) || (params[0] > range[1])) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D had bad interpolation", p);
87c1cad2e7SMatthew G. Knepley   if (Np > 1 && ((params[1] < range[2]) || (params[1] > range[3]))) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D had bad interpolation", p);
88c1cad2e7SMatthew G. Knepley   /* Put coordinates for new vertex in result[] */
89c1cad2e7SMatthew G. Knepley   ierr = EG_evaluate(obj, params, result);CHKERRQ(ierr);
90c1cad2e7SMatthew G. Knepley   for (d = 0; d < dE; ++d) gcoords[d] = result[d];
91c1cad2e7SMatthew G. Knepley   PetscFunctionReturn(0);
92c1cad2e7SMatthew G. Knepley }
93c1cad2e7SMatthew G. Knepley #endif
94c1cad2e7SMatthew G. Knepley 
95a8ededdfSMatthew G. Knepley /*@
96a8ededdfSMatthew 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.
97a8ededdfSMatthew G. Knepley 
98a8ededdfSMatthew G. Knepley   Not collective
99a8ededdfSMatthew G. Knepley 
100a8ededdfSMatthew G. Knepley   Input Parameters:
101a8ededdfSMatthew G. Knepley + dm      - The DMPlex object
102a8ededdfSMatthew G. Knepley . p       - The mesh point
103*d410b0cfSMatthew G. Knepley . dE      - The coordinate dimension
104a8ededdfSMatthew G. Knepley - mcoords - A coordinate point lying on the mesh point
105a8ededdfSMatthew G. Knepley 
106a8ededdfSMatthew G. Knepley   Output Parameter:
107a8ededdfSMatthew G. Knepley . gcoords - The closest coordinate point on the geometry model associated with 'p' to the given point
108a8ededdfSMatthew G. Knepley 
109*d410b0cfSMatthew G. Knepley   Note: Returns the original coordinates if no geometry model is found. Right now the only supported geometry model is EGADS. The coordinate dimension may be different from the coordinate dimension of the dm, for example if the transformation is extrusion.
110a8ededdfSMatthew G. Knepley 
111a8ededdfSMatthew G. Knepley   Level: intermediate
112a8ededdfSMatthew G. Knepley 
113a8ededdfSMatthew G. Knepley .seealso: DMRefine(), DMPlexCreate(), DMPlexSetRefinementUniform()
114a8ededdfSMatthew G. Knepley @*/
115*d410b0cfSMatthew G. Knepley PetscErrorCode DMPlexSnapToGeomModel(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[])
116a8ededdfSMatthew G. Knepley {
117*d410b0cfSMatthew G. Knepley   PetscInt d;
118a8ededdfSMatthew G. Knepley 
119c1cad2e7SMatthew G. Knepley   PetscFunctionBeginHot;
120a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS
121c1cad2e7SMatthew G. Knepley   {
122c1cad2e7SMatthew G. Knepley     DM_Plex       *plex = (DM_Plex *) dm->data;
123c1cad2e7SMatthew G. Knepley     DMLabel        bodyLabel, faceLabel, edgeLabel;
124c1cad2e7SMatthew G. Knepley     PetscInt       bodyID, faceID, edgeID;
125c1cad2e7SMatthew G. Knepley     PetscContainer modelObj;
126c1cad2e7SMatthew G. Knepley     ego            model;
127c1cad2e7SMatthew G. Knepley     PetscBool      islite = PETSC_FALSE;
128*d410b0cfSMatthew G. Knepley     PetscErrorCode ierr;
129c1cad2e7SMatthew G. Knepley 
130a8ededdfSMatthew G. Knepley     ierr = DMGetLabel(dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
131a8ededdfSMatthew G. Knepley     ierr = DMGetLabel(dm, "EGADS Face ID", &faceLabel);CHKERRQ(ierr);
132a8ededdfSMatthew G. Knepley     ierr = DMGetLabel(dm, "EGADS Edge ID", &edgeLabel);CHKERRQ(ierr);
133c1cad2e7SMatthew G. Knepley     if (!bodyLabel || !faceLabel || !edgeLabel || plex->ignoreModel) {
134f0fcf0adSMatthew G. Knepley       for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
135a8ededdfSMatthew G. Knepley       PetscFunctionReturn(0);
136a8ededdfSMatthew G. Knepley     }
137c1cad2e7SMatthew G. Knepley     ierr = PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
138c1cad2e7SMatthew G. Knepley     if (!modelObj) {
139c1cad2e7SMatthew G. Knepley       ierr = PetscObjectQuery((PetscObject) dm, "EGADSLite Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
140c1cad2e7SMatthew G. Knepley       islite = PETSC_TRUE;
141c1cad2e7SMatthew G. Knepley     }
142c1cad2e7SMatthew G. Knepley     if (!modelObj) {
143c1cad2e7SMatthew G. Knepley       for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
144c1cad2e7SMatthew G. Knepley       PetscFunctionReturn(0);
145c1cad2e7SMatthew G. Knepley     }
146c1cad2e7SMatthew G. Knepley     ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr);
147a8ededdfSMatthew G. Knepley     ierr = DMLabelGetValue(bodyLabel, p, &bodyID);CHKERRQ(ierr);
148a8ededdfSMatthew G. Knepley     ierr = DMLabelGetValue(faceLabel, p, &faceID);CHKERRQ(ierr);
149a8ededdfSMatthew G. Knepley     ierr = DMLabelGetValue(edgeLabel, p, &edgeID);CHKERRQ(ierr);
150c1cad2e7SMatthew G. Knepley     /* Allows for "Connective" Plex Edges present in models with multiple non-touching Entities */
151c1cad2e7SMatthew G. Knepley     if (bodyID < 0) {
152f0fcf0adSMatthew G. Knepley       for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
153f0fcf0adSMatthew G. Knepley       PetscFunctionReturn(0);
154a8ededdfSMatthew G. Knepley     }
155c1cad2e7SMatthew G. Knepley     if (islite) {ierr = DMPlexSnapToGeomModel_EGADSLite_Internal(dm, p, model, bodyID, faceID, edgeID, mcoords, gcoords);CHKERRQ(ierr);}
156c1cad2e7SMatthew G. Knepley     else        {ierr = DMPlexSnapToGeomModel_EGADS_Internal(dm, p, model, bodyID, faceID, edgeID, mcoords, gcoords);CHKERRQ(ierr);}
157f0fcf0adSMatthew G. Knepley   }
158a8ededdfSMatthew G. Knepley #else
159f0fcf0adSMatthew G. Knepley   for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
160a8ededdfSMatthew G. Knepley #endif
161a8ededdfSMatthew G. Knepley   PetscFunctionReturn(0);
162a8ededdfSMatthew G. Knepley }
1637bee2925SMatthew Knepley 
1647bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS)
165c1cad2e7SMatthew G. Knepley static PetscErrorCode DMPlexEGADSPrintModel_Internal(ego model)
1667bee2925SMatthew Knepley {
1677bee2925SMatthew Knepley   ego            geom, *bodies, *objs, *nobjs, *mobjs, *lobjs;
1687bee2925SMatthew Knepley   int            oclass, mtype, *senses;
1697bee2925SMatthew Knepley   int            Nb, b;
1707bee2925SMatthew Knepley   PetscErrorCode ierr;
1717bee2925SMatthew Knepley 
1727bee2925SMatthew Knepley   PetscFunctionBeginUser;
1737bee2925SMatthew Knepley   /* test bodyTopo functions */
1747bee2925SMatthew Knepley   ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
1757bee2925SMatthew Knepley   ierr = PetscPrintf(PETSC_COMM_SELF, " Number of BODIES (nbodies): %d \n", Nb);CHKERRQ(ierr);
1767bee2925SMatthew Knepley 
1777bee2925SMatthew Knepley   for (b = 0; b < Nb; ++b) {
1787bee2925SMatthew Knepley     ego body = bodies[b];
1797bee2925SMatthew Knepley     int id, Nsh, Nf, Nl, l, Ne, e, Nv, v;
1807bee2925SMatthew Knepley 
1817bee2925SMatthew Knepley     /* Output Basic Model Topology */
1827bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, SHELL, &Nsh, &objs);CHKERRQ(ierr);
1837bee2925SMatthew Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of SHELLS: %d \n", Nsh);CHKERRQ(ierr);
1847bee2925SMatthew Knepley     EG_free(objs);
1857bee2925SMatthew Knepley 
1867bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, FACE,  &Nf, &objs);CHKERRQ(ierr);
1877bee2925SMatthew Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of FACES: %d \n", Nf);CHKERRQ(ierr);
1887bee2925SMatthew Knepley     EG_free(objs);
1897bee2925SMatthew Knepley 
1907bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, LOOP,  &Nl, &lobjs);CHKERRQ(ierr);
1917bee2925SMatthew Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of LOOPS: %d \n", Nl);CHKERRQ(ierr);
1927bee2925SMatthew Knepley 
1937bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, EDGE,  &Ne, &objs);CHKERRQ(ierr);
1947bee2925SMatthew Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of EDGES: %d \n", Ne);CHKERRQ(ierr);
1957bee2925SMatthew Knepley     EG_free(objs);
1967bee2925SMatthew Knepley 
1977bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, NODE,  &Nv, &objs);CHKERRQ(ierr);
1987bee2925SMatthew Knepley     ierr = PetscPrintf(PETSC_COMM_SELF, "   Number of NODES: %d \n", Nv);CHKERRQ(ierr);
1997bee2925SMatthew Knepley     EG_free(objs);
2007bee2925SMatthew Knepley 
2017bee2925SMatthew Knepley     for (l = 0; l < Nl; ++l) {
2027bee2925SMatthew Knepley       ego loop = lobjs[l];
2037bee2925SMatthew Knepley 
2047bee2925SMatthew Knepley       id   = EG_indexBodyTopo(body, loop);
2057bee2925SMatthew Knepley       ierr = PetscPrintf(PETSC_COMM_SELF, "          LOOP ID: %d\n", id);CHKERRQ(ierr);
2067bee2925SMatthew Knepley 
2077bee2925SMatthew Knepley       /* Get EDGE info which associated with the current LOOP */
2087bee2925SMatthew Knepley       ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
2097bee2925SMatthew Knepley 
2107bee2925SMatthew Knepley       for (e = 0; e < Ne; ++e) {
2117bee2925SMatthew Knepley         ego    edge      = objs[e];
2127bee2925SMatthew Knepley         double range[4]  = {0., 0., 0., 0.};
2137bee2925SMatthew Knepley         double point[3]  = {0., 0., 0.};
214266cfabeSMatthew G. Knepley         double params[3] = {0., 0., 0.};
215266cfabeSMatthew G. Knepley         double result[18];
2167bee2925SMatthew Knepley         int    peri;
2177bee2925SMatthew Knepley 
2187bee2925SMatthew Knepley         id   = EG_indexBodyTopo(body, edge);CHKERRQ(ierr);
219266cfabeSMatthew G. Knepley         ierr = PetscPrintf(PETSC_COMM_SELF, "            EDGE ID: %d (%d)\n", id, e);CHKERRQ(ierr);
2207bee2925SMatthew Knepley 
221266cfabeSMatthew G. Knepley         ierr = EG_getRange(edge, range, &peri);CHKERRQ(ierr);
2221e1ea65dSPierre Jolivet         ierr = PetscPrintf(PETSC_COMM_SELF, "  Range = %lf, %lf, %lf, %lf \n", range[0], range[1], range[2], range[3]);CHKERRQ(ierr);
2237bee2925SMatthew Knepley 
2247bee2925SMatthew Knepley         /* Get NODE info which associated with the current EDGE */
2257bee2925SMatthew Knepley         ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);CHKERRQ(ierr);
226266cfabeSMatthew G. Knepley         if (mtype == DEGENERATE) {
227266cfabeSMatthew G. Knepley           ierr = PetscPrintf(PETSC_COMM_SELF, "  EDGE %d is DEGENERATE \n", id);CHKERRQ(ierr);
228266cfabeSMatthew G. Knepley         } else {
229266cfabeSMatthew G. Knepley           params[0] = range[0];
230266cfabeSMatthew G. Knepley           ierr = EG_evaluate(edge, params, result);CHKERRQ(ierr);
2311e1ea65dSPierre Jolivet           ierr = PetscPrintf(PETSC_COMM_SELF, "   between (%lf, %lf, %lf)", result[0], result[1], result[2]);CHKERRQ(ierr);
232266cfabeSMatthew G. Knepley           params[0] = range[1];
233266cfabeSMatthew G. Knepley           ierr = EG_evaluate(edge, params, result);CHKERRQ(ierr);
2341e1ea65dSPierre Jolivet           ierr = PetscPrintf(PETSC_COMM_SELF, " and (%lf, %lf, %lf)\n", result[0], result[1], result[2]);CHKERRQ(ierr);
235266cfabeSMatthew G. Knepley         }
2367bee2925SMatthew Knepley 
2377bee2925SMatthew Knepley         for (v = 0; v < Nv; ++v) {
2387bee2925SMatthew Knepley           ego    vertex = nobjs[v];
2397bee2925SMatthew Knepley           double limits[4];
2407bee2925SMatthew Knepley           int    dummy;
2417bee2925SMatthew Knepley 
2427bee2925SMatthew Knepley           ierr = EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses);CHKERRQ(ierr);
2437bee2925SMatthew Knepley           id   = EG_indexBodyTopo(body, vertex);
2447bee2925SMatthew Knepley           ierr = PetscPrintf(PETSC_COMM_SELF, "              NODE ID: %d \n", id);CHKERRQ(ierr);
2451e1ea65dSPierre Jolivet           ierr = PetscPrintf(PETSC_COMM_SELF, "                 (x, y, z) = (%lf, %lf, %lf) \n", limits[0], limits[1], limits[2]);CHKERRQ(ierr);
2467bee2925SMatthew Knepley 
2477bee2925SMatthew Knepley           point[0] = point[0] + limits[0];
2487bee2925SMatthew Knepley           point[1] = point[1] + limits[1];
2497bee2925SMatthew Knepley           point[2] = point[2] + limits[2];
2507bee2925SMatthew Knepley         }
2517bee2925SMatthew Knepley       }
2527bee2925SMatthew Knepley     }
2537bee2925SMatthew Knepley     EG_free(lobjs);
2547bee2925SMatthew Knepley   }
2557bee2925SMatthew Knepley   PetscFunctionReturn(0);
2567bee2925SMatthew Knepley }
2577bee2925SMatthew Knepley 
2587bee2925SMatthew Knepley static PetscErrorCode DMPlexEGADSDestroy_Private(void *context)
2597bee2925SMatthew Knepley {
2607bee2925SMatthew Knepley   if (context) EG_close((ego) context);
2617bee2925SMatthew Knepley   return 0;
2627bee2925SMatthew Knepley }
2637bee2925SMatthew Knepley 
264c1cad2e7SMatthew G. Knepley static PetscErrorCode DMPlexCreateEGADS_Internal(MPI_Comm comm, ego context, ego model, DM *newdm)
2657bee2925SMatthew Knepley {
266c1cad2e7SMatthew G. Knepley   DMLabel        bodyLabel, faceLabel, edgeLabel, vertexLabel;
2677bee2925SMatthew Knepley   PetscInt       cStart, cEnd, c;
2687bee2925SMatthew Knepley   /* EGADSLite variables */
2697bee2925SMatthew Knepley   ego            geom, *bodies, *objs, *nobjs, *mobjs, *lobjs;
2707bee2925SMatthew Knepley   int            oclass, mtype, nbodies, *senses;
2717bee2925SMatthew Knepley   int            b;
2727bee2925SMatthew Knepley   /* PETSc variables */
2737bee2925SMatthew Knepley   DM             dm;
2747bee2925SMatthew Knepley   PetscHMapI     edgeMap = NULL;
2757bee2925SMatthew 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;
2767bee2925SMatthew Knepley   PetscInt      *cells  = NULL, *cone = NULL;
2777bee2925SMatthew Knepley   PetscReal     *coords = NULL;
2787bee2925SMatthew Knepley   PetscMPIInt    rank;
2797bee2925SMatthew Knepley   PetscErrorCode ierr;
2807bee2925SMatthew Knepley 
2817bee2925SMatthew Knepley   PetscFunctionBegin;
28255b25c41SPierre Jolivet   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
283dd400576SPatrick Sanan   if (rank == 0) {
284266cfabeSMatthew G. Knepley     const PetscInt debug = 0;
285266cfabeSMatthew G. Knepley 
2867bee2925SMatthew Knepley     /* ---------------------------------------------------------------------------------------------------
2877bee2925SMatthew Knepley     Generate Petsc Plex
2887bee2925SMatthew Knepley       Get all Nodes in model, record coordinates in a correctly formatted array
2897bee2925SMatthew Knepley       Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array
2907bee2925SMatthew Knepley       We need to uniformly refine the initial geometry to guarantee a valid mesh
2917bee2925SMatthew Knepley     */
2927bee2925SMatthew Knepley 
2937bee2925SMatthew Knepley     /* Calculate cell and vertex sizes */
2947bee2925SMatthew Knepley     ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses);CHKERRQ(ierr);
2957bee2925SMatthew Knepley     ierr = PetscHMapICreate(&edgeMap);CHKERRQ(ierr);
2967bee2925SMatthew Knepley     numEdges = 0;
2977bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
2987bee2925SMatthew Knepley       ego body = bodies[b];
2997bee2925SMatthew Knepley       int id, Nl, l, Nv, v;
3007bee2925SMatthew Knepley 
3017bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
3027bee2925SMatthew Knepley       for (l = 0; l < Nl; ++l) {
3037bee2925SMatthew Knepley         ego loop = lobjs[l];
3047bee2925SMatthew Knepley         int Ner  = 0, Ne, e, Nc;
3057bee2925SMatthew Knepley 
3067bee2925SMatthew Knepley         ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
3077bee2925SMatthew Knepley         for (e = 0; e < Ne; ++e) {
3087bee2925SMatthew Knepley           ego edge = objs[e];
3097bee2925SMatthew Knepley           int Nv, id;
3107bee2925SMatthew Knepley           PetscHashIter iter;
3117bee2925SMatthew Knepley           PetscBool     found;
3127bee2925SMatthew Knepley 
3137bee2925SMatthew Knepley           ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);
3147bee2925SMatthew Knepley           if (mtype == DEGENERATE) continue;
3157bee2925SMatthew Knepley           id   = EG_indexBodyTopo(body, edge);CHKERRQ(ierr);
3167bee2925SMatthew Knepley           ierr = PetscHMapIFind(edgeMap, id-1, &iter, &found);CHKERRQ(ierr);
3177bee2925SMatthew Knepley           if (!found) {ierr = PetscHMapISet(edgeMap, id-1, numEdges++);CHKERRQ(ierr);}
3187bee2925SMatthew Knepley           ++Ner;
3197bee2925SMatthew Knepley         }
3207bee2925SMatthew Knepley         if (Ner == 2)      {Nc = 2;}
3217bee2925SMatthew Knepley         else if (Ner == 3) {Nc = 4;}
3227bee2925SMatthew Knepley         else if (Ner == 4) {Nc = 8; ++numQuads;}
3237bee2925SMatthew Knepley         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot support loop with %d edges", Ner);
3247bee2925SMatthew Knepley         numCells += Nc;
3257bee2925SMatthew Knepley         newCells += Nc-1;
3267bee2925SMatthew Knepley         maxCorners = PetscMax(Ner*2+1, maxCorners);
3277bee2925SMatthew Knepley       }
3287bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs);CHKERRQ(ierr);
3297bee2925SMatthew Knepley       for (v = 0; v < Nv; ++v) {
3307bee2925SMatthew Knepley         ego vertex = nobjs[v];
3317bee2925SMatthew Knepley 
3327bee2925SMatthew Knepley         id = EG_indexBodyTopo(body, vertex);
3337bee2925SMatthew Knepley         /* TODO: Instead of assuming contiguous ids, we could use a hash table */
3347bee2925SMatthew Knepley         numVertices = PetscMax(id, numVertices);
3357bee2925SMatthew Knepley       }
3367bee2925SMatthew Knepley       EG_free(lobjs);
3377bee2925SMatthew Knepley       EG_free(nobjs);
3387bee2925SMatthew Knepley     }
3397bee2925SMatthew Knepley     ierr = PetscHMapIGetSize(edgeMap, &numEdges);CHKERRQ(ierr);
3407bee2925SMatthew Knepley     newVertices  = numEdges + numQuads;
3417bee2925SMatthew Knepley     numVertices += newVertices;
3427bee2925SMatthew Knepley 
3437bee2925SMatthew Knepley     dim        = 2; /* Assume 3D Models :: Need to update to handle 2D Models in the future */
3447bee2925SMatthew Knepley     cdim       = 3; /* Assume 3D Models :: Need to update to handle 2D Models in the future */
3457bee2925SMatthew Knepley     numCorners = 3; /* Split cells into triangles */
3467bee2925SMatthew Knepley     ierr = PetscMalloc3(numVertices*cdim, &coords, numCells*numCorners, &cells, maxCorners, &cone);CHKERRQ(ierr);
3477bee2925SMatthew Knepley 
3487bee2925SMatthew Knepley     /* Get vertex coordinates */
3497bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
3507bee2925SMatthew Knepley       ego body = bodies[b];
3517bee2925SMatthew Knepley       int id, Nv, v;
3527bee2925SMatthew Knepley 
3537bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs);CHKERRQ(ierr);
3547bee2925SMatthew Knepley       for (v = 0; v < Nv; ++v) {
3557bee2925SMatthew Knepley         ego    vertex = nobjs[v];
3567bee2925SMatthew Knepley         double limits[4];
3577bee2925SMatthew Knepley         int    dummy;
3587bee2925SMatthew Knepley 
3597bee2925SMatthew Knepley         ierr = EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses);CHKERRQ(ierr);
3607bee2925SMatthew Knepley         id   = EG_indexBodyTopo(body, vertex);CHKERRQ(ierr);
3617bee2925SMatthew Knepley         coords[(id-1)*cdim+0] = limits[0];
3627bee2925SMatthew Knepley         coords[(id-1)*cdim+1] = limits[1];
3637bee2925SMatthew Knepley         coords[(id-1)*cdim+2] = limits[2];
3647bee2925SMatthew Knepley       }
3657bee2925SMatthew Knepley       EG_free(nobjs);
3667bee2925SMatthew Knepley     }
3677bee2925SMatthew Knepley     ierr = PetscHMapIClear(edgeMap);CHKERRQ(ierr);
3687bee2925SMatthew Knepley     fOff     = numVertices - newVertices + numEdges;
3697bee2925SMatthew Knepley     numEdges = 0;
3707bee2925SMatthew Knepley     numQuads = 0;
3717bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
3727bee2925SMatthew Knepley       ego body = bodies[b];
3737bee2925SMatthew Knepley       int Nl, l;
3747bee2925SMatthew Knepley 
3757bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
3767bee2925SMatthew Knepley       for (l = 0; l < Nl; ++l) {
3777bee2925SMatthew Knepley         ego loop = lobjs[l];
3787bee2925SMatthew Knepley         int lid, Ner = 0, Ne, e;
3797bee2925SMatthew Knepley 
3807bee2925SMatthew Knepley         lid  = EG_indexBodyTopo(body, loop);CHKERRQ(ierr);
3817bee2925SMatthew Knepley         ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
3827bee2925SMatthew Knepley         for (e = 0; e < Ne; ++e) {
3837bee2925SMatthew Knepley           ego       edge = objs[e];
3847bee2925SMatthew Knepley           int       eid, Nv;
3857bee2925SMatthew Knepley           PetscHashIter iter;
3867bee2925SMatthew Knepley           PetscBool     found;
3877bee2925SMatthew Knepley 
3887bee2925SMatthew Knepley           ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);
3897bee2925SMatthew Knepley           if (mtype == DEGENERATE) continue;
3907bee2925SMatthew Knepley           ++Ner;
3917bee2925SMatthew Knepley           eid  = EG_indexBodyTopo(body, edge);CHKERRQ(ierr);
3927bee2925SMatthew Knepley           ierr = PetscHMapIFind(edgeMap, eid-1, &iter, &found);CHKERRQ(ierr);
3937bee2925SMatthew Knepley           if (!found) {
3947bee2925SMatthew Knepley             PetscInt v = numVertices - newVertices + numEdges;
3957bee2925SMatthew Knepley             double range[4], params[3] = {0., 0., 0.}, result[18];
3967bee2925SMatthew Knepley             int    periodic[2];
3977bee2925SMatthew Knepley 
3987bee2925SMatthew Knepley             ierr = PetscHMapISet(edgeMap, eid-1, numEdges++);CHKERRQ(ierr);
3997bee2925SMatthew Knepley             ierr = EG_getRange(edge, range, periodic);CHKERRQ(ierr);
4007bee2925SMatthew Knepley             params[0] = 0.5*(range[0] + range[1]);
4017bee2925SMatthew Knepley             ierr = EG_evaluate(edge, params, result);CHKERRQ(ierr);
4027bee2925SMatthew Knepley             coords[v*cdim+0] = result[0];
4037bee2925SMatthew Knepley             coords[v*cdim+1] = result[1];
4047bee2925SMatthew Knepley             coords[v*cdim+2] = result[2];
4057bee2925SMatthew Knepley           }
4067bee2925SMatthew Knepley         }
4077bee2925SMatthew Knepley         if (Ner == 4) {
4087bee2925SMatthew Knepley           PetscInt v = fOff + numQuads++;
409266cfabeSMatthew G. Knepley           ego     *fobjs, face;
4107bee2925SMatthew Knepley           double   range[4], params[3] = {0., 0., 0.}, result[18];
411266cfabeSMatthew G. Knepley           int      Nf, fid, periodic[2];
4127bee2925SMatthew Knepley 
4137bee2925SMatthew Knepley           ierr = EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs);CHKERRQ(ierr);
414266cfabeSMatthew G. Knepley           face = fobjs[0];
415266cfabeSMatthew G. Knepley           fid  = EG_indexBodyTopo(body, face);CHKERRQ(ierr);
416266cfabeSMatthew 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);
417266cfabeSMatthew G. Knepley           ierr = EG_getRange(face, range, periodic);CHKERRQ(ierr);
4187bee2925SMatthew Knepley           params[0] = 0.5*(range[0] + range[1]);
4197bee2925SMatthew Knepley           params[1] = 0.5*(range[2] + range[3]);
420266cfabeSMatthew G. Knepley           ierr = EG_evaluate(face, params, result);CHKERRQ(ierr);
4217bee2925SMatthew Knepley           coords[v*cdim+0] = result[0];
4227bee2925SMatthew Knepley           coords[v*cdim+1] = result[1];
4237bee2925SMatthew Knepley           coords[v*cdim+2] = result[2];
4247bee2925SMatthew Knepley         }
4257bee2925SMatthew Knepley       }
4267bee2925SMatthew Knepley     }
4277bee2925SMatthew Knepley     if (numEdges + numQuads != newVertices) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of new vertices %D != %D previous count", numEdges + numQuads, newVertices);
4287bee2925SMatthew Knepley 
4297bee2925SMatthew Knepley     /* Get cell vertices by traversing loops */
4307bee2925SMatthew Knepley     numQuads = 0;
4317bee2925SMatthew Knepley     cOff     = 0;
4327bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
4337bee2925SMatthew Knepley       ego body = bodies[b];
4347bee2925SMatthew Knepley       int id, Nl, l;
4357bee2925SMatthew Knepley 
4367bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
4377bee2925SMatthew Knepley       for (l = 0; l < Nl; ++l) {
4387bee2925SMatthew Knepley         ego loop = lobjs[l];
4397bee2925SMatthew Knepley         int lid, Ner = 0, Ne, e, nc = 0, c, Nt, t;
4407bee2925SMatthew Knepley 
4417bee2925SMatthew Knepley         lid  = EG_indexBodyTopo(body, loop);CHKERRQ(ierr);
4427bee2925SMatthew Knepley         ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
4437bee2925SMatthew Knepley 
4447bee2925SMatthew Knepley         for (e = 0; e < Ne; ++e) {
4457bee2925SMatthew Knepley           ego edge = objs[e];
4467bee2925SMatthew Knepley           int points[3];
4477bee2925SMatthew Knepley           int eid, Nv, v, tmp;
4487bee2925SMatthew Knepley 
4497bee2925SMatthew Knepley           eid  = EG_indexBodyTopo(body, edge);
4507bee2925SMatthew Knepley           ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);
451266cfabeSMatthew G. Knepley           if (mtype == DEGENERATE) continue;
452266cfabeSMatthew G. Knepley           else                     ++Ner;
4537bee2925SMatthew Knepley           if (Nv != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Edge %d has %d vertices != 2", eid, Nv);
4547bee2925SMatthew Knepley 
4557bee2925SMatthew Knepley           for (v = 0; v < Nv; ++v) {
4567bee2925SMatthew Knepley             ego vertex = nobjs[v];
4577bee2925SMatthew Knepley 
4587bee2925SMatthew Knepley             id = EG_indexBodyTopo(body, vertex);
4597bee2925SMatthew Knepley             points[v*2] = id-1;
4607bee2925SMatthew Knepley           }
4617bee2925SMatthew Knepley           {
4627bee2925SMatthew Knepley             PetscInt edgeNum;
4637bee2925SMatthew Knepley 
4647bee2925SMatthew Knepley             ierr = PetscHMapIGet(edgeMap, eid-1, &edgeNum);CHKERRQ(ierr);
4657bee2925SMatthew Knepley             points[1] = numVertices - newVertices + edgeNum;
4667bee2925SMatthew Knepley           }
4677bee2925SMatthew Knepley           /* EGADS loops are not oriented, but seem to be in order, so we must piece them together */
4687bee2925SMatthew Knepley           if (!nc) {
4697bee2925SMatthew Knepley             for (v = 0; v < Nv+1; ++v) cone[nc++] = points[v];
4707bee2925SMatthew Knepley           } else {
4717bee2925SMatthew Knepley             if (cone[nc-1] == points[0])      {cone[nc++] = points[1]; if (cone[0] != points[2]) cone[nc++] = points[2];}
4727bee2925SMatthew Knepley             else if (cone[nc-1] == points[2]) {cone[nc++] = points[1]; if (cone[0] != points[0]) cone[nc++] = points[0];}
4737bee2925SMatthew 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];}
4747bee2925SMatthew 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];}
4757bee2925SMatthew Knepley             else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %d does not match its predecessor", eid);
4767bee2925SMatthew Knepley           }
4777bee2925SMatthew Knepley         }
4787bee2925SMatthew Knepley         if (nc != 2*Ner)     SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %D != %D", nc, 2*Ner);
4797bee2925SMatthew Knepley         if (Ner == 4) {cone[nc++] = numVertices - newVertices + numEdges + numQuads++;}
4807bee2925SMatthew Knepley         if (nc > maxCorners) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %D > %D max", nc, maxCorners);
4817bee2925SMatthew Knepley         /* Triangulate the loop */
4827bee2925SMatthew Knepley         switch (Ner) {
4837bee2925SMatthew Knepley           case 2: /* Bi-Segment -> 2 triangles */
4847bee2925SMatthew Knepley             Nt = 2;
4857bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[0];
4867bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[1];
4877bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[2];
4887bee2925SMatthew Knepley             ++cOff;
4897bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[0];
4907bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[2];
4917bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[3];
4927bee2925SMatthew Knepley             ++cOff;
4937bee2925SMatthew Knepley             break;
4947bee2925SMatthew Knepley           case 3: /* Triangle   -> 4 triangles */
4957bee2925SMatthew Knepley             Nt = 4;
4967bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[0];
4977bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[1];
4987bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[5];
4997bee2925SMatthew Knepley             ++cOff;
5007bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[1];
5017bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[2];
5027bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[3];
5037bee2925SMatthew Knepley             ++cOff;
5047bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[5];
5057bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[3];
5067bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[4];
5077bee2925SMatthew Knepley             ++cOff;
5087bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[1];
5097bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[3];
5107bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[5];
5117bee2925SMatthew Knepley             ++cOff;
5127bee2925SMatthew Knepley             break;
5137bee2925SMatthew Knepley           case 4: /* Quad       -> 8 triangles */
5147bee2925SMatthew Knepley             Nt = 8;
5157bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[0];
5167bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[1];
5177bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[7];
5187bee2925SMatthew Knepley             ++cOff;
5197bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[1];
5207bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[2];
5217bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[3];
5227bee2925SMatthew Knepley             ++cOff;
5237bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[3];
5247bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[4];
5257bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[5];
5267bee2925SMatthew Knepley             ++cOff;
5277bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[5];
5287bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[6];
5297bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[7];
5307bee2925SMatthew Knepley             ++cOff;
5317bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[8];
5327bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[1];
5337bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[3];
5347bee2925SMatthew Knepley             ++cOff;
5357bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[8];
5367bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[3];
5377bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[5];
5387bee2925SMatthew Knepley             ++cOff;
5397bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[8];
5407bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[5];
5417bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[7];
5427bee2925SMatthew Knepley             ++cOff;
5437bee2925SMatthew Knepley             cells[cOff*numCorners+0] = cone[8];
5447bee2925SMatthew Knepley             cells[cOff*numCorners+1] = cone[7];
5457bee2925SMatthew Knepley             cells[cOff*numCorners+2] = cone[1];
5467bee2925SMatthew Knepley             ++cOff;
5477bee2925SMatthew Knepley             break;
5487bee2925SMatthew Knepley           default: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d edges, which we do not support", lid, Ner);
5497bee2925SMatthew Knepley         }
550266cfabeSMatthew G. Knepley         if (debug) {
5517bee2925SMatthew Knepley           for (t = 0; t < Nt; ++t) {
5527bee2925SMatthew Knepley             ierr = PetscPrintf(PETSC_COMM_SELF, "  LOOP Corner NODEs Triangle %D (", t);CHKERRQ(ierr);
5537bee2925SMatthew Knepley             for (c = 0; c < numCorners; ++c) {
5547bee2925SMatthew Knepley               if (c > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, ", ");CHKERRQ(ierr);}
5557bee2925SMatthew Knepley               ierr = PetscPrintf(PETSC_COMM_SELF, "%D", cells[(cOff-Nt+t)*numCorners+c]);CHKERRQ(ierr);
5567bee2925SMatthew Knepley             }
5577bee2925SMatthew Knepley             ierr = PetscPrintf(PETSC_COMM_SELF, ")\n");CHKERRQ(ierr);
5587bee2925SMatthew Knepley           }
5597bee2925SMatthew Knepley         }
560266cfabeSMatthew G. Knepley       }
5617bee2925SMatthew Knepley       EG_free(lobjs);
5627bee2925SMatthew Knepley     }
5637bee2925SMatthew Knepley   }
5647bee2925SMatthew Knepley   if (cOff != numCells) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count of total cells %D != %D previous count", cOff, numCells);
5657bee2925SMatthew Knepley   ierr = DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numVertices, numCorners, PETSC_TRUE, cells, cdim, coords, &dm);CHKERRQ(ierr);
5667bee2925SMatthew Knepley   ierr = PetscFree3(coords, cells, cone);CHKERRQ(ierr);
567266cfabeSMatthew G. Knepley   ierr = PetscInfo2(dm, " Total Number of Unique Cells    = %D (%D)\n", numCells, newCells);CHKERRQ(ierr);
568266cfabeSMatthew G. Knepley   ierr = PetscInfo2(dm, " Total Number of Unique Vertices = %D (%D)\n", numVertices, newVertices);CHKERRQ(ierr);
5697bee2925SMatthew Knepley   /* Embed EGADS model in DM */
5707bee2925SMatthew Knepley   {
5717bee2925SMatthew Knepley     PetscContainer modelObj, contextObj;
5727bee2925SMatthew Knepley 
5737bee2925SMatthew Knepley     ierr = PetscContainerCreate(PETSC_COMM_SELF, &modelObj);CHKERRQ(ierr);
5747bee2925SMatthew Knepley     ierr = PetscContainerSetPointer(modelObj, model);CHKERRQ(ierr);
5757bee2925SMatthew Knepley     ierr = PetscObjectCompose((PetscObject) dm, "EGADS Model", (PetscObject) modelObj);CHKERRQ(ierr);
5767bee2925SMatthew Knepley     ierr = PetscContainerDestroy(&modelObj);CHKERRQ(ierr);
5777bee2925SMatthew Knepley 
5787bee2925SMatthew Knepley     ierr = PetscContainerCreate(PETSC_COMM_SELF, &contextObj);CHKERRQ(ierr);
5797bee2925SMatthew Knepley     ierr = PetscContainerSetPointer(contextObj, context);CHKERRQ(ierr);
5807bee2925SMatthew Knepley     ierr = PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private);CHKERRQ(ierr);
5817bee2925SMatthew Knepley     ierr = PetscObjectCompose((PetscObject) dm, "EGADS Context", (PetscObject) contextObj);CHKERRQ(ierr);
5827bee2925SMatthew Knepley     ierr = PetscContainerDestroy(&contextObj);CHKERRQ(ierr);
5837bee2925SMatthew Knepley   }
5847bee2925SMatthew Knepley   /* Label points */
5857bee2925SMatthew Knepley   ierr = DMCreateLabel(dm, "EGADS Body ID");CHKERRQ(ierr);
5867bee2925SMatthew Knepley   ierr = DMGetLabel(dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
5877bee2925SMatthew Knepley   ierr = DMCreateLabel(dm, "EGADS Face ID");CHKERRQ(ierr);
5887bee2925SMatthew Knepley   ierr = DMGetLabel(dm, "EGADS Face ID", &faceLabel);CHKERRQ(ierr);
5897bee2925SMatthew Knepley   ierr = DMCreateLabel(dm, "EGADS Edge ID");CHKERRQ(ierr);
5907bee2925SMatthew Knepley   ierr = DMGetLabel(dm, "EGADS Edge ID", &edgeLabel);CHKERRQ(ierr);
591c1cad2e7SMatthew G. Knepley   ierr = DMCreateLabel(dm, "EGADS Vertex ID");CHKERRQ(ierr);
592c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel);CHKERRQ(ierr);
5937bee2925SMatthew Knepley   cOff = 0;
5947bee2925SMatthew Knepley   for (b = 0; b < nbodies; ++b) {
5957bee2925SMatthew Knepley     ego body = bodies[b];
5967bee2925SMatthew Knepley     int id, Nl, l;
5977bee2925SMatthew Knepley 
5987bee2925SMatthew Knepley     ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
5997bee2925SMatthew Knepley     for (l = 0; l < Nl; ++l) {
6007bee2925SMatthew Knepley       ego  loop = lobjs[l];
6017bee2925SMatthew Knepley       ego *fobjs;
6027bee2925SMatthew Knepley       int  lid, Nf, fid, Ner = 0, Ne, e, Nt = 0, t;
6037bee2925SMatthew Knepley 
6047bee2925SMatthew Knepley       lid  = EG_indexBodyTopo(body, loop);CHKERRQ(ierr);
6057bee2925SMatthew Knepley       ierr = EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs);CHKERRQ(ierr);
6062758c9b9SBarry Smith       if (Nf > 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf);
6077bee2925SMatthew Knepley       fid  = EG_indexBodyTopo(body, fobjs[0]);CHKERRQ(ierr);
6087bee2925SMatthew Knepley       EG_free(fobjs);
6097bee2925SMatthew Knepley       ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr);
6107bee2925SMatthew Knepley       for (e = 0; e < Ne; ++e) {
6117bee2925SMatthew Knepley         ego             edge = objs[e];
6127bee2925SMatthew Knepley         int             eid, Nv, v;
6137bee2925SMatthew Knepley         PetscInt        points[3], support[2], numEdges, edgeNum;
6147bee2925SMatthew Knepley         const PetscInt *edges;
6157bee2925SMatthew Knepley 
6167bee2925SMatthew Knepley         eid  = EG_indexBodyTopo(body, edge);
6177bee2925SMatthew Knepley         ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);
6187bee2925SMatthew Knepley         if (mtype == DEGENERATE) continue;
6197bee2925SMatthew Knepley         else                     ++Ner;
6207bee2925SMatthew Knepley         for (v = 0; v < Nv; ++v) {
6217bee2925SMatthew Knepley           ego vertex = nobjs[v];
6227bee2925SMatthew Knepley 
6237bee2925SMatthew Knepley           id   = EG_indexBodyTopo(body, vertex);
6247bee2925SMatthew Knepley           ierr = DMLabelSetValue(edgeLabel, numCells + id-1, eid);CHKERRQ(ierr);
6257bee2925SMatthew Knepley           points[v*2] = numCells + id-1;
6267bee2925SMatthew Knepley         }
6277bee2925SMatthew Knepley         ierr = PetscHMapIGet(edgeMap, eid-1, &edgeNum);CHKERRQ(ierr);
6287bee2925SMatthew Knepley         points[1] = numCells + numVertices - newVertices + edgeNum;
6297bee2925SMatthew Knepley 
6307bee2925SMatthew Knepley         ierr = DMLabelSetValue(edgeLabel, points[1], eid);CHKERRQ(ierr);
6317bee2925SMatthew Knepley         support[0] = points[0];
6327bee2925SMatthew Knepley         support[1] = points[1];
6337bee2925SMatthew Knepley         ierr = DMPlexGetJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr);
6347bee2925SMatthew 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);
6357bee2925SMatthew Knepley         ierr = DMLabelSetValue(edgeLabel, edges[0], eid);CHKERRQ(ierr);
6367bee2925SMatthew Knepley         ierr = DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr);
6377bee2925SMatthew Knepley         support[0] = points[1];
6387bee2925SMatthew Knepley         support[1] = points[2];
6397bee2925SMatthew Knepley         ierr = DMPlexGetJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr);
6407bee2925SMatthew 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);
6417bee2925SMatthew Knepley         ierr = DMLabelSetValue(edgeLabel, edges[0], eid);CHKERRQ(ierr);
6427bee2925SMatthew Knepley         ierr = DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr);
6437bee2925SMatthew Knepley       }
6447bee2925SMatthew Knepley       switch (Ner) {
6457bee2925SMatthew Knepley         case 2: Nt = 2;break;
6467bee2925SMatthew Knepley         case 3: Nt = 4;break;
6477bee2925SMatthew Knepley         case 4: Nt = 8;break;
6487bee2925SMatthew Knepley         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Loop with %d edges is unsupported", Ner);
6497bee2925SMatthew Knepley       }
6507bee2925SMatthew Knepley       for (t = 0; t < Nt; ++t) {
6517bee2925SMatthew Knepley         ierr = DMLabelSetValue(bodyLabel, cOff+t, b);CHKERRQ(ierr);
6527bee2925SMatthew Knepley         ierr = DMLabelSetValue(faceLabel, cOff+t, fid);CHKERRQ(ierr);
6537bee2925SMatthew Knepley       }
6547bee2925SMatthew Knepley       cOff += Nt;
6557bee2925SMatthew Knepley     }
6567bee2925SMatthew Knepley     EG_free(lobjs);
6577bee2925SMatthew Knepley   }
6587bee2925SMatthew Knepley   ierr = PetscHMapIDestroy(&edgeMap);CHKERRQ(ierr);
6597bee2925SMatthew Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
6607bee2925SMatthew Knepley   for (c = cStart; c < cEnd; ++c) {
6617bee2925SMatthew Knepley     PetscInt *closure = NULL;
6627bee2925SMatthew Knepley     PetscInt  clSize, cl, bval, fval;
6637bee2925SMatthew Knepley 
6647bee2925SMatthew Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
6657bee2925SMatthew Knepley     ierr = DMLabelGetValue(bodyLabel, c, &bval);CHKERRQ(ierr);
6667bee2925SMatthew Knepley     ierr = DMLabelGetValue(faceLabel, c, &fval);CHKERRQ(ierr);
6677bee2925SMatthew Knepley     for (cl = 0; cl < clSize*2; cl += 2) {
6687bee2925SMatthew Knepley       ierr = DMLabelSetValue(bodyLabel, closure[cl], bval);CHKERRQ(ierr);
6697bee2925SMatthew Knepley       ierr = DMLabelSetValue(faceLabel, closure[cl], fval);CHKERRQ(ierr);
6707bee2925SMatthew Knepley     }
6717bee2925SMatthew Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr);
6727bee2925SMatthew Knepley   }
6737bee2925SMatthew Knepley   *newdm = dm;
6747bee2925SMatthew Knepley   PetscFunctionReturn(0);
6757bee2925SMatthew Knepley }
676c1cad2e7SMatthew G. Knepley 
677c1cad2e7SMatthew G. Knepley static PetscErrorCode DMPlexCreateEGADS(MPI_Comm comm, ego context, ego model, DM *newdm)
678c1cad2e7SMatthew G. Knepley {
679c1cad2e7SMatthew G. Knepley   DMLabel         bodyLabel, faceLabel, edgeLabel, vertexLabel;
680c1cad2e7SMatthew G. Knepley   // EGADS/EGADSLite variables
681c1cad2e7SMatthew G. Knepley   ego             geom, *bodies, *mobjs, *fobjs, *lobjs, *eobjs, *nobjs;
682c1cad2e7SMatthew G. Knepley   ego             topRef, prev, next;
683c1cad2e7SMatthew G. Knepley   int             oclass, mtype, nbodies, *senses, *lSenses, *eSenses;
684c1cad2e7SMatthew G. Knepley   int             b;
685c1cad2e7SMatthew G. Knepley   // PETSc variables
686c1cad2e7SMatthew G. Knepley   DM              dm;
687c1cad2e7SMatthew G. Knepley   PetscHMapI      edgeMap = NULL, bodyIndexMap = NULL, bodyVertexMap = NULL, bodyEdgeMap = NULL, bodyFaceMap = NULL, bodyEdgeGlobalMap = NULL;
688c1cad2e7SMatthew G. Knepley   PetscInt        dim = -1, cdim = -1, numCorners = 0, numVertices = 0, numEdges = 0, numFaces = 0, numCells = 0, edgeCntr = 0;
689c1cad2e7SMatthew G. Knepley   PetscInt        cellCntr = 0, numPoints = 0;
690c1cad2e7SMatthew G. Knepley   PetscInt        *cells  = NULL;
691c1cad2e7SMatthew G. Knepley   const PetscInt  *cone = NULL;
692c1cad2e7SMatthew G. Knepley   PetscReal       *coords = NULL;
693c1cad2e7SMatthew G. Knepley   PetscMPIInt      rank;
694c1cad2e7SMatthew G. Knepley   PetscErrorCode   ierr;
695c1cad2e7SMatthew G. Knepley 
696c1cad2e7SMatthew G. Knepley   PetscFunctionBeginUser;
697c1cad2e7SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
698c1cad2e7SMatthew G. Knepley   if (!rank) {
699c1cad2e7SMatthew G. Knepley     // ---------------------------------------------------------------------------------------------------
700c1cad2e7SMatthew G. Knepley     // Generate Petsc Plex
701c1cad2e7SMatthew G. Knepley     //  Get all Nodes in model, record coordinates in a correctly formatted array
702c1cad2e7SMatthew G. Knepley     //  Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array
703c1cad2e7SMatthew G. Knepley     //  We need to uniformly refine the initial geometry to guarantee a valid mesh
704c1cad2e7SMatthew G. Knepley 
705c1cad2e7SMatthew G. Knepley   // Caluculate cell and vertex sizes
706c1cad2e7SMatthew G. Knepley   ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses);CHKERRQ(ierr);
707c1cad2e7SMatthew G. Knepley 
708c1cad2e7SMatthew G. Knepley     ierr = PetscHMapICreate(&edgeMap);CHKERRQ(ierr);
709c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&bodyIndexMap);CHKERRQ(ierr);
710c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&bodyVertexMap);CHKERRQ(ierr);
711c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&bodyEdgeMap);CHKERRQ(ierr);
712c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&bodyEdgeGlobalMap);CHKERRQ(ierr);
713c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&bodyFaceMap);CHKERRQ(ierr);
714c1cad2e7SMatthew G. Knepley 
715c1cad2e7SMatthew G. Knepley   for (b = 0; b < nbodies; ++b) {
716c1cad2e7SMatthew G. Knepley       ego             body = bodies[b];
717c1cad2e7SMatthew G. Knepley     int             Nf, Ne, Nv;
718c1cad2e7SMatthew G. Knepley     PetscHashIter   BIiter, BViter, BEiter, BEGiter, BFiter, EMiter;
719c1cad2e7SMatthew G. Knepley     PetscBool       BIfound, BVfound, BEfound, BEGfound, BFfound, EMfound;
720c1cad2e7SMatthew G. Knepley 
721c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyIndexMap, b, &BIiter, &BIfound);CHKERRQ(ierr);
722c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound);CHKERRQ(ierr);
723c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound);CHKERRQ(ierr);
724c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound);CHKERRQ(ierr);
725c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound);CHKERRQ(ierr);
726c1cad2e7SMatthew G. Knepley 
727c1cad2e7SMatthew G. Knepley     if (!BIfound)  {ierr = PetscHMapISet(bodyIndexMap, b, numFaces + numEdges + numVertices);CHKERRQ(ierr);}
728c1cad2e7SMatthew G. Knepley     if (!BVfound)  {ierr = PetscHMapISet(bodyVertexMap, b, numVertices);CHKERRQ(ierr);}
729c1cad2e7SMatthew G. Knepley     if (!BEfound)  {ierr = PetscHMapISet(bodyEdgeMap, b, numEdges);CHKERRQ(ierr);}
730c1cad2e7SMatthew G. Knepley     if (!BEGfound) {ierr = PetscHMapISet(bodyEdgeGlobalMap, b, edgeCntr);CHKERRQ(ierr);}
731c1cad2e7SMatthew G. Knepley     if (!BFfound)  {ierr = PetscHMapISet(bodyFaceMap, b, numFaces);CHKERRQ(ierr);}
732c1cad2e7SMatthew G. Knepley 
733c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs);CHKERRQ(ierr);
734c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs);CHKERRQ(ierr);
735c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs);CHKERRQ(ierr);
736c1cad2e7SMatthew G. Knepley     EG_free(fobjs);
737c1cad2e7SMatthew G. Knepley     EG_free(eobjs);
738c1cad2e7SMatthew G. Knepley     EG_free(nobjs);
739c1cad2e7SMatthew G. Knepley 
740c1cad2e7SMatthew G. Knepley     // Remove DEGENERATE EDGES from Edge count
741c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs);CHKERRQ(ierr);
742c1cad2e7SMatthew G. Knepley     int Netemp = 0;
743c1cad2e7SMatthew G. Knepley     for (int e = 0; e < Ne; ++e) {
744c1cad2e7SMatthew G. Knepley       ego     edge = eobjs[e];
745c1cad2e7SMatthew G. Knepley       int     eid;
746c1cad2e7SMatthew G. Knepley 
747c1cad2e7SMatthew G. Knepley       ierr = EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next);CHKERRQ(ierr);
748c1cad2e7SMatthew G. Knepley       eid = EG_indexBodyTopo(body, edge);CHKERRQ(ierr);
749c1cad2e7SMatthew G. Knepley 
750c1cad2e7SMatthew G. Knepley       ierr = PetscHMapIFind(edgeMap, edgeCntr + eid - 1, &EMiter, &EMfound);CHKERRQ(ierr);
751c1cad2e7SMatthew G. Knepley       if (mtype == DEGENERATE) {
752c1cad2e7SMatthew G. Knepley         if (!EMfound) {ierr = PetscHMapISet(edgeMap, edgeCntr + eid - 1, -1);CHKERRQ(ierr);}
753c1cad2e7SMatthew G. Knepley       }
754c1cad2e7SMatthew G. Knepley       else {
755c1cad2e7SMatthew G. Knepley       ++Netemp;
756c1cad2e7SMatthew G. Knepley         if (!EMfound) {ierr = PetscHMapISet(edgeMap, edgeCntr + eid - 1, Netemp);CHKERRQ(ierr);}
757c1cad2e7SMatthew G. Knepley       }
758c1cad2e7SMatthew G. Knepley     }
759c1cad2e7SMatthew G. Knepley     EG_free(eobjs);
760c1cad2e7SMatthew G. Knepley 
761c1cad2e7SMatthew G. Knepley     // Determine Number of Cells
762c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs);CHKERRQ(ierr);
763c1cad2e7SMatthew G. Knepley     for (int f = 0; f < Nf; ++f) {
764c1cad2e7SMatthew G. Knepley         ego     face = fobjs[f];
765c1cad2e7SMatthew G. Knepley     int     edgeTemp = 0;
766c1cad2e7SMatthew G. Knepley 
767c1cad2e7SMatthew G. Knepley       ierr = EG_getBodyTopos(body, face, EDGE, &Ne, &eobjs);CHKERRQ(ierr);
768c1cad2e7SMatthew G. Knepley       for (int e = 0; e < Ne; ++e) {
769c1cad2e7SMatthew G. Knepley         ego     edge = eobjs[e];
770c1cad2e7SMatthew G. Knepley 
771c1cad2e7SMatthew G. Knepley         ierr = EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next);CHKERRQ(ierr);
772c1cad2e7SMatthew G. Knepley         if (mtype != DEGENERATE) {++edgeTemp;}
773c1cad2e7SMatthew G. Knepley       }
774c1cad2e7SMatthew G. Knepley       numCells += (2 * edgeTemp);
775c1cad2e7SMatthew G. Knepley       EG_free(eobjs);
776c1cad2e7SMatthew G. Knepley     }
777c1cad2e7SMatthew G. Knepley     EG_free(fobjs);
778c1cad2e7SMatthew G. Knepley 
779c1cad2e7SMatthew G. Knepley     numFaces    += Nf;
780c1cad2e7SMatthew G. Knepley     numEdges    += Netemp;
781c1cad2e7SMatthew G. Knepley     numVertices += Nv;
782c1cad2e7SMatthew G. Knepley     edgeCntr    += Ne;
783c1cad2e7SMatthew G. Knepley   }
784c1cad2e7SMatthew G. Knepley 
785c1cad2e7SMatthew G. Knepley   // Set up basic DMPlex parameters
786c1cad2e7SMatthew G. Knepley   dim        = 2;    // Assumes 3D Models :: Need to handle 2D modles in the future
787c1cad2e7SMatthew G. Knepley   cdim       = 3;     // Assumes 3D Models :: Need to update to handle 2D modles in future
788c1cad2e7SMatthew G. Knepley   numCorners = 3;     // Split Faces into triangles
789c1cad2e7SMatthew G. Knepley     numPoints  = numVertices + numEdges + numFaces;   // total number of coordinate points
790c1cad2e7SMatthew G. Knepley 
791c1cad2e7SMatthew G. Knepley   ierr = PetscMalloc2(numPoints*cdim, &coords, numCells*numCorners, &cells);CHKERRQ(ierr);
792c1cad2e7SMatthew G. Knepley 
793c1cad2e7SMatthew G. Knepley   // Get Vertex Coordinates and Set up Cells
794c1cad2e7SMatthew G. Knepley   for (b = 0; b < nbodies; ++b) {
795c1cad2e7SMatthew G. Knepley     ego             body = bodies[b];
796c1cad2e7SMatthew G. Knepley     int             Nf, Ne, Nv;
797c1cad2e7SMatthew G. Knepley     PetscInt        bodyVertexIndexStart, bodyEdgeIndexStart, bodyEdgeGlobalIndexStart, bodyFaceIndexStart;
798c1cad2e7SMatthew G. Knepley     PetscHashIter   BViter, BEiter, BEGiter, BFiter, EMiter;
799c1cad2e7SMatthew G. Knepley     PetscBool       BVfound, BEfound, BEGfound, BFfound, EMfound;
800c1cad2e7SMatthew G. Knepley 
801c1cad2e7SMatthew G. Knepley     // Vertices on Current Body
802c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs);CHKERRQ(ierr);
803c1cad2e7SMatthew G. Knepley 
804c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound);CHKERRQ(ierr);
805c1cad2e7SMatthew G. Knepley     if (!BVfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyVertexMap", b);
806c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIGet(bodyVertexMap, b, &bodyVertexIndexStart);CHKERRQ(ierr);
807c1cad2e7SMatthew G. Knepley 
808c1cad2e7SMatthew G. Knepley     for (int v = 0; v < Nv; ++v) {
809c1cad2e7SMatthew G. Knepley       ego    vertex = nobjs[v];
810c1cad2e7SMatthew G. Knepley     double limits[4];
811c1cad2e7SMatthew G. Knepley     int    id, dummy;
812c1cad2e7SMatthew G. Knepley 
813c1cad2e7SMatthew G. Knepley     ierr = EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses);CHKERRQ(ierr);
814c1cad2e7SMatthew G. Knepley     id = EG_indexBodyTopo(body, vertex);CHKERRQ(ierr);
815c1cad2e7SMatthew G. Knepley 
816c1cad2e7SMatthew G. Knepley     coords[(bodyVertexIndexStart + id - 1)*cdim + 0] = limits[0];
817c1cad2e7SMatthew G. Knepley     coords[(bodyVertexIndexStart + id - 1)*cdim + 1] = limits[1];
818c1cad2e7SMatthew G. Knepley     coords[(bodyVertexIndexStart + id - 1)*cdim + 2] = limits[2];
819c1cad2e7SMatthew G. Knepley     }
820c1cad2e7SMatthew G. Knepley     EG_free(nobjs);
821c1cad2e7SMatthew G. Knepley 
822c1cad2e7SMatthew G. Knepley     // Edge Midpoint Vertices on Current Body
823c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs);CHKERRQ(ierr);
824c1cad2e7SMatthew G. Knepley 
825c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound);CHKERRQ(ierr);
826c1cad2e7SMatthew G. Knepley     if (!BEfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeMap", b);
827c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIGet(bodyEdgeMap, b, &bodyEdgeIndexStart);CHKERRQ(ierr);
828c1cad2e7SMatthew G. Knepley 
829c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound);CHKERRQ(ierr);
830c1cad2e7SMatthew G. Knepley     if (!BEGfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeGlobalMap", b);
831c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIGet(bodyEdgeGlobalMap, b, &bodyEdgeGlobalIndexStart);CHKERRQ(ierr);
832c1cad2e7SMatthew G. Knepley 
833c1cad2e7SMatthew G. Knepley     for (int e = 0; e < Ne; ++e) {
834c1cad2e7SMatthew G. Knepley       ego          edge = eobjs[e];
835c1cad2e7SMatthew G. Knepley     double       range[2], avgt[1], cntrPnt[9];
836c1cad2e7SMatthew G. Knepley     int          eid, eOffset;
837c1cad2e7SMatthew G. Knepley     int          periodic;
838c1cad2e7SMatthew G. Knepley 
839c1cad2e7SMatthew G. Knepley     ierr = EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next);CHKERRQ(ierr);
840c1cad2e7SMatthew G. Knepley     if (mtype == DEGENERATE) {continue;}
841c1cad2e7SMatthew G. Knepley 
842c1cad2e7SMatthew G. Knepley     eid = EG_indexBodyTopo(body, edge);CHKERRQ(ierr);
843c1cad2e7SMatthew G. Knepley 
844c1cad2e7SMatthew G. Knepley     // get relative offset from globalEdgeID Vector
845c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound);CHKERRQ(ierr);
846c1cad2e7SMatthew G. Knepley       if (!EMfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %d not found in edgeMap", bodyEdgeGlobalIndexStart + eid - 1);
847c1cad2e7SMatthew G. Knepley       ierr = PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset);CHKERRQ(ierr);
848c1cad2e7SMatthew G. Knepley 
849c1cad2e7SMatthew G. Knepley     ierr = EG_getRange(edge, range, &periodic);CHKERRQ(ierr);
850c1cad2e7SMatthew G. Knepley     avgt[0] = (range[0] + range[1]) /  2.;
851c1cad2e7SMatthew G. Knepley 
852c1cad2e7SMatthew G. Knepley     ierr = EG_evaluate(edge, avgt, cntrPnt);CHKERRQ(ierr);
853c1cad2e7SMatthew G. Knepley     coords[(numVertices + bodyEdgeIndexStart + eOffset - 1)*cdim + 0] = cntrPnt[0];
854c1cad2e7SMatthew G. Knepley         coords[(numVertices + bodyEdgeIndexStart + eOffset - 1)*cdim + 1] = cntrPnt[1];
855c1cad2e7SMatthew G. Knepley     coords[(numVertices + bodyEdgeIndexStart + eOffset - 1)*cdim + 2] = cntrPnt[2];
856c1cad2e7SMatthew G. Knepley     }
857c1cad2e7SMatthew G. Knepley     EG_free(eobjs);
858c1cad2e7SMatthew G. Knepley 
859c1cad2e7SMatthew G. Knepley     // Face Midpoint Vertices on Current Body
860c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs);CHKERRQ(ierr);
861c1cad2e7SMatthew G. Knepley 
862c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound);CHKERRQ(ierr);
863c1cad2e7SMatthew G. Knepley     if (!BFfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyFaceMap", b);
864c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIGet(bodyFaceMap, b, &bodyFaceIndexStart);CHKERRQ(ierr);
865c1cad2e7SMatthew G. Knepley 
866c1cad2e7SMatthew G. Knepley     for (int f = 0; f < Nf; ++f) {
867c1cad2e7SMatthew G. Knepley     ego       face = fobjs[f];
868c1cad2e7SMatthew G. Knepley     double    range[4], avgUV[2], cntrPnt[18];
869c1cad2e7SMatthew G. Knepley     int       peri, id;
870c1cad2e7SMatthew G. Knepley 
871c1cad2e7SMatthew G. Knepley     id = EG_indexBodyTopo(body, face);
872c1cad2e7SMatthew G. Knepley     ierr = EG_getRange(face, range, &peri);CHKERRQ(ierr);
873c1cad2e7SMatthew G. Knepley 
874c1cad2e7SMatthew G. Knepley     avgUV[0] = (range[0] + range[1]) / 2.;
875c1cad2e7SMatthew G. Knepley     avgUV[1] = (range[2] + range[3]) / 2.;
876c1cad2e7SMatthew G. Knepley     ierr = EG_evaluate(face, avgUV, cntrPnt);CHKERRQ(ierr);
877c1cad2e7SMatthew G. Knepley 
878c1cad2e7SMatthew G. Knepley     coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1)*cdim + 0] = cntrPnt[0];
879c1cad2e7SMatthew G. Knepley     coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1)*cdim + 1] = cntrPnt[1];
880c1cad2e7SMatthew G. Knepley     coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1)*cdim + 2] = cntrPnt[2];
881c1cad2e7SMatthew G. Knepley     }
882c1cad2e7SMatthew G. Knepley     EG_free(fobjs);
883c1cad2e7SMatthew G. Knepley 
884c1cad2e7SMatthew G. Knepley     // Define Cells :: Note - This could be incorporated in the Face Midpoint Vertices Loop but was kept separate for clarity
885c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs);CHKERRQ(ierr);
886c1cad2e7SMatthew G. Knepley     for (int f = 0; f < Nf; ++f) {
887c1cad2e7SMatthew G. Knepley     ego      face = fobjs[f];
888c1cad2e7SMatthew G. Knepley     int      fID, midFaceID, midPntID, startID, endID, Nl;
889c1cad2e7SMatthew G. Knepley 
890c1cad2e7SMatthew G. Knepley     fID = EG_indexBodyTopo(body, face);CHKERRQ(ierr);
891c1cad2e7SMatthew G. Knepley     midFaceID = numVertices + numEdges + bodyFaceIndexStart + fID - 1;
892c1cad2e7SMatthew G. Knepley     // Must Traverse Loop to ensure we have all necessary information like the sense (+/- 1) of the edges.
893c1cad2e7SMatthew G. Knepley     // TODO :: Only handles single loop faces (No holes). The choices for handling multiloop faces are:
894c1cad2e7SMatthew G. Knepley     //            1) Use the DMPlexCreateEGADSFromFile() with the -dm_plex_egads_with_tess = 1 option.
895c1cad2e7SMatthew G. Knepley     //               This will use a default EGADS tessellation as an initial surface mesh.
896c1cad2e7SMatthew G. Knepley     //            2) Create the initial surface mesh via a 2D mesher :: Currently not availble (?future?)
897c1cad2e7SMatthew G. Knepley     //               May I suggest the XXXX as a starting point?
898c1cad2e7SMatthew G. Knepley 
899c1cad2e7SMatthew G. Knepley     ierr = EG_getTopology(face, &geom, &oclass, &mtype, NULL, &Nl, &lobjs, &lSenses);CHKERRQ(ierr);
900c1cad2e7SMatthew G. Knepley 
901c1cad2e7SMatthew G. Knepley       if (Nl > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Face has %d Loops. Can only handle Faces with 1 Loop. Please use --dm_plex_egads_with_tess = 1 Option", Nl);
902c1cad2e7SMatthew G. Knepley     for (int l = 0; l < Nl; ++l) {
903c1cad2e7SMatthew G. Knepley           ego      loop = lobjs[l];
904c1cad2e7SMatthew G. Knepley 
905c1cad2e7SMatthew G. Knepley           ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses);CHKERRQ(ierr);
906c1cad2e7SMatthew G. Knepley       for (int e = 0; e < Ne; ++e) {
907c1cad2e7SMatthew G. Knepley         ego     edge = eobjs[e];
908c1cad2e7SMatthew G. Knepley         int     eid, eOffset;
909c1cad2e7SMatthew G. Knepley 
910c1cad2e7SMatthew G. Knepley         ierr = EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next);CHKERRQ(ierr);
911c1cad2e7SMatthew G. Knepley       eid = EG_indexBodyTopo(body, edge);
912c1cad2e7SMatthew G. Knepley         if (mtype == DEGENERATE) { continue; }
913c1cad2e7SMatthew G. Knepley 
914c1cad2e7SMatthew G. Knepley         // get relative offset from globalEdgeID Vector
915c1cad2e7SMatthew G. Knepley         ierr = PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound);CHKERRQ(ierr);
916c1cad2e7SMatthew G. Knepley           if (!EMfound) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %d of Body %d not found in edgeMap. Global Edge ID :: %d", eid, b, bodyEdgeGlobalIndexStart + eid - 1);
917c1cad2e7SMatthew G. Knepley           ierr = PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset);CHKERRQ(ierr);
918c1cad2e7SMatthew G. Knepley 
919c1cad2e7SMatthew G. Knepley       midPntID = numVertices + bodyEdgeIndexStart + eOffset - 1;
920c1cad2e7SMatthew G. Knepley 
921c1cad2e7SMatthew G. Knepley         ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);CHKERRQ(ierr);
922c1cad2e7SMatthew G. Knepley 
923c1cad2e7SMatthew G. Knepley         if (eSenses[e] > 0) { startID = EG_indexBodyTopo(body, nobjs[0]); endID = EG_indexBodyTopo(body, nobjs[1]); }
924c1cad2e7SMatthew G. Knepley         else { startID = EG_indexBodyTopo(body, nobjs[1]); endID = EG_indexBodyTopo(body, nobjs[0]); }
925c1cad2e7SMatthew G. Knepley 
926c1cad2e7SMatthew G. Knepley       // Define 2 Cells per Edge with correct orientation
927c1cad2e7SMatthew G. Knepley       cells[cellCntr*numCorners + 0] = midFaceID;
928c1cad2e7SMatthew G. Knepley       cells[cellCntr*numCorners + 1] = bodyVertexIndexStart + startID - 1;
929c1cad2e7SMatthew G. Knepley       cells[cellCntr*numCorners + 2] = midPntID;
930c1cad2e7SMatthew G. Knepley 
931c1cad2e7SMatthew G. Knepley       cells[cellCntr*numCorners + 3] = midFaceID;
932c1cad2e7SMatthew G. Knepley       cells[cellCntr*numCorners + 4] = midPntID;
933c1cad2e7SMatthew G. Knepley       cells[cellCntr*numCorners + 5] = bodyVertexIndexStart + endID - 1;
934c1cad2e7SMatthew G. Knepley 
935c1cad2e7SMatthew G. Knepley       cellCntr = cellCntr + 2;
936c1cad2e7SMatthew G. Knepley       }
937c1cad2e7SMatthew G. Knepley     }
938c1cad2e7SMatthew G. Knepley     }
939c1cad2e7SMatthew G. Knepley     EG_free(fobjs);
940c1cad2e7SMatthew G. Knepley   }
941c1cad2e7SMatthew G. Knepley   }
942c1cad2e7SMatthew G. Knepley 
943c1cad2e7SMatthew G. Knepley   // Generate DMPlex
944c1cad2e7SMatthew G. Knepley   ierr = DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numPoints, numCorners, PETSC_TRUE, cells, cdim, coords, &dm);CHKERRQ(ierr);
945c1cad2e7SMatthew G. Knepley   ierr = PetscFree2(coords, cells);CHKERRQ(ierr);
946c1cad2e7SMatthew G. Knepley   ierr = PetscInfo1(dm, " Total Number of Unique Cells    = %D \n", numCells);CHKERRQ(ierr);
947c1cad2e7SMatthew G. Knepley   ierr = PetscInfo1(dm, " Total Number of Unique Vertices = %D \n", numVertices);CHKERRQ(ierr);
948c1cad2e7SMatthew G. Knepley 
949c1cad2e7SMatthew G. Knepley   // Embed EGADS model in DM
950c1cad2e7SMatthew G. Knepley   {
951c1cad2e7SMatthew G. Knepley     PetscContainer modelObj, contextObj;
952c1cad2e7SMatthew G. Knepley 
953c1cad2e7SMatthew G. Knepley     ierr = PetscContainerCreate(PETSC_COMM_SELF, &modelObj);CHKERRQ(ierr);
954c1cad2e7SMatthew G. Knepley     ierr = PetscContainerSetPointer(modelObj, model);CHKERRQ(ierr);
955c1cad2e7SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dm, "EGADS Model", (PetscObject) modelObj);CHKERRQ(ierr);
956c1cad2e7SMatthew G. Knepley     ierr = PetscContainerDestroy(&modelObj);CHKERRQ(ierr);
957c1cad2e7SMatthew G. Knepley 
958c1cad2e7SMatthew G. Knepley     ierr = PetscContainerCreate(PETSC_COMM_SELF, &contextObj);CHKERRQ(ierr);
959c1cad2e7SMatthew G. Knepley     ierr = PetscContainerSetPointer(contextObj, context);CHKERRQ(ierr);
960c1cad2e7SMatthew G. Knepley     ierr = PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private);CHKERRQ(ierr);
961c1cad2e7SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dm, "EGADS Context", (PetscObject) contextObj);CHKERRQ(ierr);
962c1cad2e7SMatthew G. Knepley     ierr = PetscContainerDestroy(&contextObj);CHKERRQ(ierr);
963c1cad2e7SMatthew G. Knepley   }
964c1cad2e7SMatthew G. Knepley   // Label points
965c1cad2e7SMatthew G. Knepley   PetscInt   nStart, nEnd;
966c1cad2e7SMatthew G. Knepley 
967c1cad2e7SMatthew G. Knepley   ierr = DMCreateLabel(dm, "EGADS Body ID");CHKERRQ(ierr);
968c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
969c1cad2e7SMatthew G. Knepley   ierr = DMCreateLabel(dm, "EGADS Face ID");CHKERRQ(ierr);
970c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Face ID", &faceLabel);CHKERRQ(ierr);
971c1cad2e7SMatthew G. Knepley   ierr = DMCreateLabel(dm, "EGADS Edge ID");CHKERRQ(ierr);
972c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Edge ID", &edgeLabel);CHKERRQ(ierr);
973c1cad2e7SMatthew G. Knepley   ierr = DMCreateLabel(dm, "EGADS Vertex ID");CHKERRQ(ierr);
974c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel);CHKERRQ(ierr);
975c1cad2e7SMatthew G. Knepley 
976c1cad2e7SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 2, &nStart, &nEnd);CHKERRQ(ierr);
977c1cad2e7SMatthew G. Knepley 
978c1cad2e7SMatthew G. Knepley   cellCntr = 0;
979c1cad2e7SMatthew G. Knepley   for (b = 0; b < nbodies; ++b) {
980c1cad2e7SMatthew G. Knepley     ego             body = bodies[b];
981c1cad2e7SMatthew G. Knepley   int             Nv, Ne, Nf;
982c1cad2e7SMatthew G. Knepley   PetscInt        bodyVertexIndexStart, bodyEdgeIndexStart, bodyEdgeGlobalIndexStart, bodyFaceIndexStart;
983c1cad2e7SMatthew G. Knepley   PetscHashIter   BViter, BEiter, BEGiter, BFiter, EMiter;
984c1cad2e7SMatthew G. Knepley   PetscBool       BVfound, BEfound, BEGfound, BFfound, EMfound;
985c1cad2e7SMatthew G. Knepley 
986c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound);CHKERRQ(ierr);
987c1cad2e7SMatthew G. Knepley   if (!BVfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyVertexMap", b);
988c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIGet(bodyVertexMap, b, &bodyVertexIndexStart);CHKERRQ(ierr);
989c1cad2e7SMatthew G. Knepley 
990c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound);CHKERRQ(ierr);
991c1cad2e7SMatthew G. Knepley   if (!BEfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeMap", b);
992c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIGet(bodyEdgeMap, b, &bodyEdgeIndexStart);CHKERRQ(ierr);
993c1cad2e7SMatthew G. Knepley 
994c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound);CHKERRQ(ierr);
995c1cad2e7SMatthew G. Knepley   if (!BFfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyFaceMap", b);
996c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIGet(bodyFaceMap, b, &bodyFaceIndexStart);CHKERRQ(ierr);
997c1cad2e7SMatthew G. Knepley 
998c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound);CHKERRQ(ierr);
999c1cad2e7SMatthew G. Knepley     if (!BEGfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeGlobalMap", b);
1000c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIGet(bodyEdgeGlobalMap, b, &bodyEdgeGlobalIndexStart);CHKERRQ(ierr);
1001c1cad2e7SMatthew G. Knepley 
1002c1cad2e7SMatthew G. Knepley   ierr = EG_getBodyTopos(body, NULL, FACE,  &Nf, &fobjs);CHKERRQ(ierr);
1003c1cad2e7SMatthew G. Knepley   for (int f = 0; f < Nf; ++f) {
1004c1cad2e7SMatthew G. Knepley     ego   face = fobjs[f];
1005c1cad2e7SMatthew G. Knepley       int   fID, Nl;
1006c1cad2e7SMatthew G. Knepley 
1007c1cad2e7SMatthew G. Knepley     fID  = EG_indexBodyTopo(body, face);CHKERRQ(ierr);
1008c1cad2e7SMatthew G. Knepley 
1009c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, face, LOOP, &Nl, &lobjs);CHKERRQ(ierr);
1010c1cad2e7SMatthew G. Knepley     for (int l = 0; l < Nl; ++l) {
1011c1cad2e7SMatthew G. Knepley         ego  loop = lobjs[l];
1012c1cad2e7SMatthew G. Knepley     int  lid;
1013c1cad2e7SMatthew G. Knepley 
1014c1cad2e7SMatthew G. Knepley     lid  = EG_indexBodyTopo(body, loop);CHKERRQ(ierr);
1015c1cad2e7SMatthew G. Knepley       if (Nl > 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf);
1016c1cad2e7SMatthew G. Knepley 
1017c1cad2e7SMatthew G. Knepley     ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses);CHKERRQ(ierr);
1018c1cad2e7SMatthew G. Knepley     for (int e = 0; e < Ne; ++e) {
1019c1cad2e7SMatthew G. Knepley       ego     edge = eobjs[e];
1020c1cad2e7SMatthew G. Knepley       int     eid, eOffset;
1021c1cad2e7SMatthew G. Knepley 
1022c1cad2e7SMatthew G. Knepley       // Skip DEGENERATE Edges
1023c1cad2e7SMatthew G. Knepley       ierr = EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next);CHKERRQ(ierr);
1024c1cad2e7SMatthew G. Knepley       if (mtype == DEGENERATE) {continue;}
1025c1cad2e7SMatthew G. Knepley       eid = EG_indexBodyTopo(body, edge);CHKERRQ(ierr);
1026c1cad2e7SMatthew G. Knepley 
1027c1cad2e7SMatthew G. Knepley       // get relative offset from globalEdgeID Vector
1028c1cad2e7SMatthew G. Knepley       ierr = PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound);CHKERRQ(ierr);
1029c1cad2e7SMatthew G. Knepley       if (!EMfound) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %d of Body %d not found in edgeMap. Global Edge ID :: %d", eid, b, bodyEdgeGlobalIndexStart + eid - 1);
1030c1cad2e7SMatthew G. Knepley       ierr = PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset);CHKERRQ(ierr);
1031c1cad2e7SMatthew G. Knepley 
1032c1cad2e7SMatthew G. Knepley       ierr = EG_getBodyTopos(body, edge, NODE, &Nv, &nobjs);CHKERRQ(ierr);
1033c1cad2e7SMatthew G. Knepley       for (int v = 0; v < Nv; ++v){
1034c1cad2e7SMatthew G. Knepley         ego vertex = nobjs[v];
1035c1cad2e7SMatthew G. Knepley         int vID;
1036c1cad2e7SMatthew G. Knepley 
1037c1cad2e7SMatthew G. Knepley         vID = EG_indexBodyTopo(body, vertex);CHKERRQ(ierr);
1038c1cad2e7SMatthew G. Knepley         ierr = DMLabelSetValue(bodyLabel, nStart + bodyVertexIndexStart + vID - 1, b);CHKERRQ(ierr);
1039c1cad2e7SMatthew G. Knepley         ierr = DMLabelSetValue(vertexLabel, nStart + bodyVertexIndexStart + vID - 1, vID);CHKERRQ(ierr);
1040c1cad2e7SMatthew G. Knepley       }
1041c1cad2e7SMatthew G. Knepley       EG_free(nobjs);
1042c1cad2e7SMatthew G. Knepley 
1043c1cad2e7SMatthew G. Knepley       ierr = DMLabelSetValue(bodyLabel, nStart + numVertices + bodyEdgeIndexStart + eOffset - 1, b);CHKERRQ(ierr);
1044c1cad2e7SMatthew G. Knepley       ierr = DMLabelSetValue(edgeLabel, nStart + numVertices + bodyEdgeIndexStart + eOffset - 1, eid);CHKERRQ(ierr);
1045c1cad2e7SMatthew G. Knepley 
1046c1cad2e7SMatthew G. Knepley       // Define Cell faces
1047c1cad2e7SMatthew G. Knepley       for (int jj = 0; jj < 2; ++jj){
1048c1cad2e7SMatthew G. Knepley         ierr = DMLabelSetValue(bodyLabel, cellCntr, b);CHKERRQ(ierr);
1049c1cad2e7SMatthew G. Knepley         ierr = DMLabelSetValue(faceLabel, cellCntr, fID);CHKERRQ(ierr);
1050c1cad2e7SMatthew G. Knepley         ierr = DMPlexGetCone(dm, cellCntr, &cone);CHKERRQ(ierr);
1051c1cad2e7SMatthew G. Knepley 
1052c1cad2e7SMatthew G. Knepley         ierr = DMLabelSetValue(bodyLabel, cone[0], b);CHKERRQ(ierr);
1053c1cad2e7SMatthew G. Knepley         ierr = DMLabelSetValue(faceLabel, cone[0], fID);CHKERRQ(ierr);
1054c1cad2e7SMatthew G. Knepley 
1055c1cad2e7SMatthew G. Knepley         ierr = DMLabelSetValue(bodyLabel, cone[1], b);CHKERRQ(ierr);
1056c1cad2e7SMatthew G. Knepley         ierr = DMLabelSetValue(edgeLabel, cone[1], eid);CHKERRQ(ierr);
1057c1cad2e7SMatthew G. Knepley 
1058c1cad2e7SMatthew G. Knepley        ierr = DMLabelSetValue(bodyLabel, cone[2], b);CHKERRQ(ierr);
1059c1cad2e7SMatthew G. Knepley        ierr = DMLabelSetValue(faceLabel, cone[2], fID);CHKERRQ(ierr);
1060c1cad2e7SMatthew G. Knepley 
1061c1cad2e7SMatthew G. Knepley        cellCntr = cellCntr + 1;
1062c1cad2e7SMatthew G. Knepley       }
1063c1cad2e7SMatthew G. Knepley     }
1064c1cad2e7SMatthew G. Knepley     }
1065c1cad2e7SMatthew G. Knepley     EG_free(lobjs);
1066c1cad2e7SMatthew G. Knepley 
1067c1cad2e7SMatthew G. Knepley     ierr = DMLabelSetValue(bodyLabel, nStart + numVertices + numEdges + bodyFaceIndexStart + fID - 1, b);CHKERRQ(ierr);
1068c1cad2e7SMatthew G. Knepley     ierr = DMLabelSetValue(faceLabel, nStart + numVertices + numEdges + bodyFaceIndexStart + fID - 1, fID);CHKERRQ(ierr);
1069c1cad2e7SMatthew G. Knepley   }
1070c1cad2e7SMatthew G. Knepley   EG_free(fobjs);
1071c1cad2e7SMatthew G. Knepley   }
1072c1cad2e7SMatthew G. Knepley 
1073c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIDestroy(&edgeMap);CHKERRQ(ierr);
1074c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIDestroy(&bodyIndexMap);CHKERRQ(ierr);
1075c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIDestroy(&bodyVertexMap);CHKERRQ(ierr);
1076c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIDestroy(&bodyEdgeMap);CHKERRQ(ierr);
1077c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIDestroy(&bodyEdgeGlobalMap);CHKERRQ(ierr);
1078c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIDestroy(&bodyFaceMap);CHKERRQ(ierr);
1079c1cad2e7SMatthew G. Knepley 
1080c1cad2e7SMatthew G. Knepley   *newdm = dm;
1081c1cad2e7SMatthew G. Knepley   PetscFunctionReturn(0);
1082c1cad2e7SMatthew G. Knepley }
1083c1cad2e7SMatthew G. Knepley 
1084c1cad2e7SMatthew G. Knepley static PetscErrorCode DMPlexCreateEGADS_Tess_Internal(MPI_Comm comm, ego context, ego model, DM *newdm)
1085c1cad2e7SMatthew G. Knepley {
1086c1cad2e7SMatthew G. Knepley   DMLabel              bodyLabel, faceLabel, edgeLabel, vertexLabel;
1087c1cad2e7SMatthew G. Knepley   /* EGADSLite variables */
1088c1cad2e7SMatthew G. Knepley   ego                  geom, *bodies, *fobjs;
1089c1cad2e7SMatthew G. Knepley   int                  b, oclass, mtype, nbodies, *senses;
1090c1cad2e7SMatthew G. Knepley   int                  totalNumTris = 0, totalNumPoints = 0;
1091c1cad2e7SMatthew G. Knepley   double               boundBox[6] = {0., 0., 0., 0., 0., 0.}, tessSize;
1092c1cad2e7SMatthew G. Knepley   /* PETSc variables */
1093c1cad2e7SMatthew G. Knepley   DM                   dm;
1094c1cad2e7SMatthew G. Knepley   PetscHMapI           pointIndexStartMap = NULL, triIndexStartMap = NULL, pTypeLabelMap = NULL, pIndexLabelMap = NULL;
1095c1cad2e7SMatthew G. Knepley   PetscHMapI           pBodyIndexLabelMap = NULL, triFaceIDLabelMap = NULL, triBodyIDLabelMap = NULL;
1096c1cad2e7SMatthew G. Knepley   PetscInt             dim = -1, cdim = -1, numCorners = 0, counter = 0;
1097c1cad2e7SMatthew G. Knepley   PetscInt            *cells  = NULL;
1098c1cad2e7SMatthew G. Knepley   const PetscInt      *cone = NULL;
1099c1cad2e7SMatthew G. Knepley   PetscReal           *coords = NULL;
1100c1cad2e7SMatthew G. Knepley   PetscMPIInt          rank;
1101c1cad2e7SMatthew G. Knepley   PetscErrorCode       ierr;
1102c1cad2e7SMatthew G. Knepley 
1103c1cad2e7SMatthew G. Knepley   PetscFunctionBeginUser;
1104c1cad2e7SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
1105c1cad2e7SMatthew G. Knepley   if (!rank) {
1106c1cad2e7SMatthew G. Knepley     // ---------------------------------------------------------------------------------------------------
1107c1cad2e7SMatthew G. Knepley     // Generate Petsc Plex from EGADSlite created Tessellation of geometry
1108c1cad2e7SMatthew G. Knepley     // ---------------------------------------------------------------------------------------------------
1109c1cad2e7SMatthew G. Knepley 
1110c1cad2e7SMatthew G. Knepley   // Caluculate cell and vertex sizes
1111c1cad2e7SMatthew G. Knepley   ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses);CHKERRQ(ierr);
1112c1cad2e7SMatthew G. Knepley 
1113c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&pointIndexStartMap);CHKERRQ(ierr);
1114c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&triIndexStartMap);CHKERRQ(ierr);
1115c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&pTypeLabelMap);CHKERRQ(ierr);
1116c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&pIndexLabelMap);CHKERRQ(ierr);
1117c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&pBodyIndexLabelMap);CHKERRQ(ierr);
1118c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&triFaceIDLabelMap);CHKERRQ(ierr);
1119c1cad2e7SMatthew G. Knepley   ierr = PetscHMapICreate(&triBodyIDLabelMap);CHKERRQ(ierr);
1120c1cad2e7SMatthew G. Knepley 
1121c1cad2e7SMatthew G. Knepley   /* Create Tessellation of Bodies */
1122c1cad2e7SMatthew G. Knepley   ego tessArray[nbodies];
1123c1cad2e7SMatthew G. Knepley 
1124c1cad2e7SMatthew G. Knepley   for (b = 0; b < nbodies; ++b) {
1125c1cad2e7SMatthew G. Knepley     ego             body = bodies[b];
1126c1cad2e7SMatthew G. Knepley     double          params[3] = {0.0, 0.0, 0.0};    // Parameters for Tessellation
1127c1cad2e7SMatthew G. Knepley     int             Nf, bodyNumPoints = 0, bodyNumTris = 0;
1128c1cad2e7SMatthew G. Knepley     PetscHashIter   PISiter, TISiter;
1129c1cad2e7SMatthew G. Knepley     PetscBool       PISfound, TISfound;
1130c1cad2e7SMatthew G. Knepley 
1131c1cad2e7SMatthew G. Knepley     /* Store Start Index for each Body's Point and Tris */
1132c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(pointIndexStartMap, b, &PISiter, &PISfound);CHKERRQ(ierr);
1133c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(triIndexStartMap, b, &TISiter, &TISfound);CHKERRQ(ierr);
1134c1cad2e7SMatthew G. Knepley 
1135c1cad2e7SMatthew G. Knepley     if (!PISfound)  {ierr = PetscHMapISet(pointIndexStartMap, b, totalNumPoints);CHKERRQ(ierr);}
1136c1cad2e7SMatthew G. Knepley     if (!TISfound)  {ierr = PetscHMapISet(triIndexStartMap, b, totalNumTris);CHKERRQ(ierr);}
1137c1cad2e7SMatthew G. Knepley 
1138c1cad2e7SMatthew G. Knepley     /* Calculate Tessellation parameters based on Bounding Box */
1139c1cad2e7SMatthew G. Knepley     /* Get Bounding Box Dimensions of the BODY */
1140c1cad2e7SMatthew G. Knepley     ierr = EG_getBoundingBox(body, boundBox);
1141c1cad2e7SMatthew G. Knepley     tessSize = boundBox[3] - boundBox[0];
1142c1cad2e7SMatthew G. Knepley     if (tessSize < boundBox[4] - boundBox[1]) tessSize = boundBox[4] - boundBox[1];
1143c1cad2e7SMatthew G. Knepley     if (tessSize < boundBox[5] - boundBox[2]) tessSize = boundBox[5] - boundBox[2];
1144c1cad2e7SMatthew G. Knepley 
1145c1cad2e7SMatthew G. Knepley     // TODO :: May want to give users tessellation parameter options //
1146c1cad2e7SMatthew G. Knepley     params[0] = 0.0250 * tessSize;
1147c1cad2e7SMatthew G. Knepley     params[1] = 0.0075 * tessSize;
1148c1cad2e7SMatthew G. Knepley     params[2] = 15.0;
1149c1cad2e7SMatthew G. Knepley 
1150c1cad2e7SMatthew G. Knepley     ierr = EG_makeTessBody(body, params, &tessArray[b]);CHKERRQ(ierr);
1151c1cad2e7SMatthew G. Knepley 
1152c1cad2e7SMatthew G. Knepley     ierr = EG_getBodyTopos(body, NULL, FACE,  &Nf, &fobjs);CHKERRQ(ierr);
1153c1cad2e7SMatthew G. Knepley 
1154c1cad2e7SMatthew G. Knepley     for (int f = 0; f < Nf; ++f) {
1155c1cad2e7SMatthew G. Knepley       ego             face = fobjs[f];
1156c1cad2e7SMatthew G. Knepley     int             len, fID, ntris;
1157c1cad2e7SMatthew G. Knepley     const int      *ptype, *pindex, *ptris, *ptric;
1158c1cad2e7SMatthew G. Knepley     const double   *pxyz, *puv;
1159c1cad2e7SMatthew G. Knepley 
1160c1cad2e7SMatthew G. Knepley     // Get Face ID //
1161c1cad2e7SMatthew G. Knepley     fID = EG_indexBodyTopo(body, face);
1162c1cad2e7SMatthew G. Knepley 
1163c1cad2e7SMatthew G. Knepley     // Checkout the Surface Tessellation //
1164c1cad2e7SMatthew G. Knepley     ierr = EG_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric);CHKERRQ(ierr);
1165c1cad2e7SMatthew G. Knepley 
1166c1cad2e7SMatthew G. Knepley     // Determine total number of triangle cells in the tessellation //
1167c1cad2e7SMatthew G. Knepley     bodyNumTris += (int) ntris;
1168c1cad2e7SMatthew G. Knepley 
1169c1cad2e7SMatthew G. Knepley     // Check out the point index and coordinate //
1170c1cad2e7SMatthew G. Knepley     for (int p = 0; p < len; ++p) {
1171c1cad2e7SMatthew G. Knepley         int global;
1172c1cad2e7SMatthew G. Knepley 
1173c1cad2e7SMatthew G. Knepley       ierr = EG_localToGlobal(tessArray[b], fID, p+1, &global);
1174c1cad2e7SMatthew G. Knepley 
1175c1cad2e7SMatthew G. Knepley       // Determine the total number of points in the tessellation //
1176c1cad2e7SMatthew G. Knepley         bodyNumPoints = PetscMax(bodyNumPoints, global);
1177c1cad2e7SMatthew G. Knepley     }
1178c1cad2e7SMatthew G. Knepley     }
1179c1cad2e7SMatthew G. Knepley     EG_free(fobjs);
1180c1cad2e7SMatthew G. Knepley 
1181c1cad2e7SMatthew G. Knepley     totalNumPoints += bodyNumPoints;
1182c1cad2e7SMatthew G. Knepley     totalNumTris += bodyNumTris;
1183c1cad2e7SMatthew G. Knepley     }
1184c1cad2e7SMatthew G. Knepley   //}  - Original End of (!rank)
1185c1cad2e7SMatthew G. Knepley 
1186c1cad2e7SMatthew G. Knepley   dim = 2;
1187c1cad2e7SMatthew G. Knepley   cdim = 3;
1188c1cad2e7SMatthew G. Knepley   numCorners = 3;
1189c1cad2e7SMatthew G. Knepley   //PetscInt counter = 0;
1190c1cad2e7SMatthew G. Knepley 
1191c1cad2e7SMatthew G. Knepley   /* NEED TO DEFINE MATRICES/VECTORS TO STORE GEOM REFERENCE DATA   */
1192c1cad2e7SMatthew G. Knepley   /* Fill in below and use to define DMLabels after DMPlex creation */
1193c1cad2e7SMatthew G. Knepley   ierr = PetscMalloc2(totalNumPoints*cdim, &coords, totalNumTris*numCorners, &cells);CHKERRQ(ierr);
1194c1cad2e7SMatthew G. Knepley 
1195c1cad2e7SMatthew G. Knepley   for (b = 0; b < nbodies; ++b) {
1196c1cad2e7SMatthew G. Knepley   ego             body = bodies[b];
1197c1cad2e7SMatthew G. Knepley   int             Nf;
1198c1cad2e7SMatthew G. Knepley   PetscInt        pointIndexStart;
1199c1cad2e7SMatthew G. Knepley   PetscHashIter   PISiter;
1200c1cad2e7SMatthew G. Knepley   PetscBool       PISfound;
1201c1cad2e7SMatthew G. Knepley 
1202c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIFind(pointIndexStartMap, b, &PISiter, &PISfound);CHKERRQ(ierr);
1203c1cad2e7SMatthew G. Knepley   if (!PISfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in pointIndexStartMap", b);
1204c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIGet(pointIndexStartMap, b, &pointIndexStart);CHKERRQ(ierr);
1205c1cad2e7SMatthew G. Knepley 
1206c1cad2e7SMatthew G. Knepley   ierr = EG_getBodyTopos(body, NULL, FACE,  &Nf, &fobjs);CHKERRQ(ierr);
1207c1cad2e7SMatthew G. Knepley 
1208c1cad2e7SMatthew G. Knepley   for (int f = 0; f < Nf; ++f) {
1209c1cad2e7SMatthew G. Knepley     /* Get Face Object */
1210c1cad2e7SMatthew G. Knepley     ego              face = fobjs[f];
1211c1cad2e7SMatthew G. Knepley     int              len, fID, ntris;
1212c1cad2e7SMatthew G. Knepley     const int       *ptype, *pindex, *ptris, *ptric;
1213c1cad2e7SMatthew G. Knepley     const double    *pxyz, *puv;
1214c1cad2e7SMatthew G. Knepley 
1215c1cad2e7SMatthew G. Knepley     /* Get Face ID */
1216c1cad2e7SMatthew G. Knepley     fID = EG_indexBodyTopo(body, face);
1217c1cad2e7SMatthew G. Knepley 
1218c1cad2e7SMatthew G. Knepley     /* Checkout the Surface Tessellation */
1219c1cad2e7SMatthew G. Knepley     ierr = EG_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric);CHKERRQ(ierr);
1220c1cad2e7SMatthew G. Knepley 
1221c1cad2e7SMatthew G. Knepley     /* Check out the point index and coordinate */
1222c1cad2e7SMatthew G. Knepley     for (int p = 0; p < len; ++p) {
1223c1cad2e7SMatthew G. Knepley     int              global;
1224c1cad2e7SMatthew G. Knepley     PetscHashIter    PTLiter, PILiter, PBLiter;
1225c1cad2e7SMatthew G. Knepley     PetscBool        PTLfound, PILfound, PBLfound;
1226c1cad2e7SMatthew G. Knepley 
1227c1cad2e7SMatthew G. Knepley     ierr = EG_localToGlobal(tessArray[b], fID, p+1, &global);
1228c1cad2e7SMatthew G. Knepley 
1229c1cad2e7SMatthew G. Knepley     /* Set the coordinates array for DAG */
1230c1cad2e7SMatthew G. Knepley     coords[((global-1+pointIndexStart)*3) + 0] = pxyz[(p*3)+0];
1231c1cad2e7SMatthew G. Knepley     coords[((global-1+pointIndexStart)*3) + 1] = pxyz[(p*3)+1];
1232c1cad2e7SMatthew G. Knepley     coords[((global-1+pointIndexStart)*3) + 2] = pxyz[(p*3)+2];
1233c1cad2e7SMatthew G. Knepley 
1234c1cad2e7SMatthew G. Knepley     /* Store Geometry Label Information for DMLabel assignment later */
1235c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(pTypeLabelMap, global-1+pointIndexStart, &PTLiter, &PTLfound);CHKERRQ(ierr);
1236c1cad2e7SMatthew G. Knepley         ierr = PetscHMapIFind(pIndexLabelMap, global-1+pointIndexStart, &PILiter, &PILfound);CHKERRQ(ierr);
1237c1cad2e7SMatthew G. Knepley         ierr = PetscHMapIFind(pBodyIndexLabelMap, global-1+pointIndexStart, &PBLiter, &PBLfound);CHKERRQ(ierr);
1238c1cad2e7SMatthew G. Knepley 
1239c1cad2e7SMatthew G. Knepley         if (!PTLfound)  {ierr = PetscHMapISet(pTypeLabelMap, global-1+pointIndexStart, ptype[p]);CHKERRQ(ierr);}
1240c1cad2e7SMatthew G. Knepley         if (!PILfound)  {ierr = PetscHMapISet(pIndexLabelMap, global-1+pointIndexStart, pindex[p]);CHKERRQ(ierr);}
1241c1cad2e7SMatthew G. Knepley         if (!PBLfound)  {ierr = PetscHMapISet(pBodyIndexLabelMap, global-1+pointIndexStart, b);CHKERRQ(ierr);}
1242c1cad2e7SMatthew G. Knepley 
1243c1cad2e7SMatthew G. Knepley     if (ptype[p] < 0) { ierr = PetscHMapISet(pIndexLabelMap, global-1+pointIndexStart, fID);CHKERRQ(ierr);}
1244c1cad2e7SMatthew G. Knepley     }
1245c1cad2e7SMatthew G. Knepley 
1246c1cad2e7SMatthew G. Knepley     for (int t = 0; t < (int) ntris; ++t){
1247c1cad2e7SMatthew G. Knepley     int             global, globalA, globalB;
1248c1cad2e7SMatthew G. Knepley     PetscHashIter   TFLiter, TBLiter;
1249c1cad2e7SMatthew G. Knepley       PetscBool       TFLfound, TBLfound;
1250c1cad2e7SMatthew G. Knepley 
1251c1cad2e7SMatthew G. Knepley     ierr = EG_localToGlobal(tessArray[b], fID, ptris[(t*3) + 0], &global);
1252c1cad2e7SMatthew G. Knepley     cells[(counter*3) +0] = global-1+pointIndexStart;
1253c1cad2e7SMatthew G. Knepley 
1254c1cad2e7SMatthew G. Knepley     ierr = EG_localToGlobal(tessArray[b], fID, ptris[(t*3) + 1], &globalA);
1255c1cad2e7SMatthew G. Knepley     cells[(counter*3) +1] = globalA-1+pointIndexStart;
1256c1cad2e7SMatthew G. Knepley 
1257c1cad2e7SMatthew G. Knepley     ierr = EG_localToGlobal(tessArray[b], fID, ptris[(t*3) + 2], &globalB);
1258c1cad2e7SMatthew G. Knepley     cells[(counter*3) +2] = globalB-1+pointIndexStart;
1259c1cad2e7SMatthew G. Knepley 
1260c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(triFaceIDLabelMap, counter, &TFLiter, &TFLfound);CHKERRQ(ierr);
1261c1cad2e7SMatthew G. Knepley         ierr = PetscHMapIFind(triBodyIDLabelMap, counter, &TBLiter, &TBLfound);CHKERRQ(ierr);
1262c1cad2e7SMatthew G. Knepley 
1263c1cad2e7SMatthew G. Knepley     if (!TFLfound)  {ierr = PetscHMapISet(triFaceIDLabelMap, counter, fID);CHKERRQ(ierr);}
1264c1cad2e7SMatthew G. Knepley         if (!TBLfound)  {ierr = PetscHMapISet(triBodyIDLabelMap, counter, b);CHKERRQ(ierr);}
1265c1cad2e7SMatthew G. Knepley 
1266c1cad2e7SMatthew G. Knepley     counter += 1;
1267c1cad2e7SMatthew G. Knepley     }
1268c1cad2e7SMatthew G. Knepley   }
1269c1cad2e7SMatthew G. Knepley   EG_free(fobjs);
1270c1cad2e7SMatthew G. Knepley   }
1271c1cad2e7SMatthew G. Knepley }
1272c1cad2e7SMatthew G. Knepley 
1273c1cad2e7SMatthew G. Knepley   //Build DMPlex
1274c1cad2e7SMatthew G. Knepley   ierr = DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, totalNumTris, totalNumPoints, numCorners, PETSC_TRUE, cells, cdim, coords, &dm);CHKERRQ(ierr);
1275c1cad2e7SMatthew G. Knepley   ierr = PetscFree2(coords, cells);CHKERRQ(ierr);
1276c1cad2e7SMatthew G. Knepley 
1277c1cad2e7SMatthew G. Knepley   // Embed EGADS model in DM
1278c1cad2e7SMatthew G. Knepley   {
1279c1cad2e7SMatthew G. Knepley     PetscContainer modelObj, contextObj;
1280c1cad2e7SMatthew G. Knepley 
1281c1cad2e7SMatthew G. Knepley     ierr = PetscContainerCreate(PETSC_COMM_SELF, &modelObj);CHKERRQ(ierr);
1282c1cad2e7SMatthew G. Knepley     ierr = PetscContainerSetPointer(modelObj, model);CHKERRQ(ierr);
1283c1cad2e7SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dm, "EGADS Model", (PetscObject) modelObj);CHKERRQ(ierr);
1284c1cad2e7SMatthew G. Knepley     ierr = PetscContainerDestroy(&modelObj);CHKERRQ(ierr);
1285c1cad2e7SMatthew G. Knepley 
1286c1cad2e7SMatthew G. Knepley     ierr = PetscContainerCreate(PETSC_COMM_SELF, &contextObj);CHKERRQ(ierr);
1287c1cad2e7SMatthew G. Knepley     ierr = PetscContainerSetPointer(contextObj, context);CHKERRQ(ierr);
1288c1cad2e7SMatthew G. Knepley     ierr = PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private);CHKERRQ(ierr);
1289c1cad2e7SMatthew G. Knepley     ierr = PetscObjectCompose((PetscObject) dm, "EGADS Context", (PetscObject) contextObj);CHKERRQ(ierr);
1290c1cad2e7SMatthew G. Knepley     ierr = PetscContainerDestroy(&contextObj);CHKERRQ(ierr);
1291c1cad2e7SMatthew G. Knepley   }
1292c1cad2e7SMatthew G. Knepley 
1293c1cad2e7SMatthew G. Knepley   // Label Points
1294c1cad2e7SMatthew G. Knepley   ierr = DMCreateLabel(dm, "EGADS Body ID");CHKERRQ(ierr);
1295c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
1296c1cad2e7SMatthew G. Knepley   ierr = DMCreateLabel(dm, "EGADS Face ID");CHKERRQ(ierr);
1297c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Face ID", &faceLabel);CHKERRQ(ierr);
1298c1cad2e7SMatthew G. Knepley   ierr = DMCreateLabel(dm, "EGADS Edge ID");CHKERRQ(ierr);
1299c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Edge ID", &edgeLabel);CHKERRQ(ierr);
1300c1cad2e7SMatthew G. Knepley   ierr = DMCreateLabel(dm, "EGADS Vertex ID");CHKERRQ(ierr);
1301c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel);CHKERRQ(ierr);
1302c1cad2e7SMatthew G. Knepley 
1303c1cad2e7SMatthew G. Knepley    /* Get Number of DAG Nodes at each level */
1304c1cad2e7SMatthew G. Knepley   int   fStart, fEnd, eStart, eEnd, nStart, nEnd;
1305c1cad2e7SMatthew G. Knepley 
1306c1cad2e7SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &fStart, &fEnd);CHKERRQ(ierr);
1307c1cad2e7SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
1308c1cad2e7SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 2, &nStart, &nEnd);CHKERRQ(ierr);
1309c1cad2e7SMatthew G. Knepley 
1310c1cad2e7SMatthew G. Knepley   /* Set DMLabels for NODES */
1311c1cad2e7SMatthew G. Knepley   for (int n = nStart; n < nEnd; ++n) {
1312c1cad2e7SMatthew G. Knepley     int             pTypeVal, pIndexVal, pBodyVal;
1313c1cad2e7SMatthew G. Knepley     PetscHashIter   PTLiter, PILiter, PBLiter;
1314c1cad2e7SMatthew G. Knepley     PetscBool       PTLfound, PILfound, PBLfound;
1315c1cad2e7SMatthew G. Knepley 
1316c1cad2e7SMatthew G. Knepley     //Converted to Hash Tables
1317c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(pTypeLabelMap, n - nStart, &PTLiter, &PTLfound);CHKERRQ(ierr);
1318c1cad2e7SMatthew G. Knepley     if (!PTLfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pTypeLabelMap", n);
1319c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIGet(pTypeLabelMap, n - nStart, &pTypeVal);CHKERRQ(ierr);
1320c1cad2e7SMatthew G. Knepley 
1321c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(pIndexLabelMap, n - nStart, &PILiter, &PILfound);CHKERRQ(ierr);
1322c1cad2e7SMatthew G. Knepley     if (!PILfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pIndexLabelMap", n);
1323c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIGet(pIndexLabelMap, n - nStart, &pIndexVal);CHKERRQ(ierr);
1324c1cad2e7SMatthew G. Knepley 
1325c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIFind(pBodyIndexLabelMap, n - nStart, &PBLiter, &PBLfound);CHKERRQ(ierr);
1326c1cad2e7SMatthew G. Knepley     if (!PBLfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pBodyLabelMap", n);
1327c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIGet(pBodyIndexLabelMap, n - nStart, &pBodyVal);CHKERRQ(ierr);
1328c1cad2e7SMatthew G. Knepley 
1329c1cad2e7SMatthew G. Knepley     ierr = DMLabelSetValue(bodyLabel, n, pBodyVal);CHKERRQ(ierr);
1330c1cad2e7SMatthew G. Knepley     if (pTypeVal == 0) {ierr = DMLabelSetValue(vertexLabel, n, pIndexVal);CHKERRQ(ierr);}
1331c1cad2e7SMatthew G. Knepley     if (pTypeVal >  0) {ierr = DMLabelSetValue(edgeLabel, n, pIndexVal);CHKERRQ(ierr);}
1332c1cad2e7SMatthew G. Knepley     if (pTypeVal <  0) {ierr = DMLabelSetValue(faceLabel, n, pIndexVal);CHKERRQ(ierr);}
1333c1cad2e7SMatthew G. Knepley   }
1334c1cad2e7SMatthew G. Knepley 
1335c1cad2e7SMatthew G. Knepley   /* Set DMLabels for Edges - Based on the DMLabels of the EDGE's NODES */
1336c1cad2e7SMatthew G. Knepley   for (int e = eStart; e < eEnd; ++e) {
1337c1cad2e7SMatthew G. Knepley   int    bodyID_0, vertexID_0, vertexID_1, edgeID_0, edgeID_1, faceID_0, faceID_1;
1338c1cad2e7SMatthew G. Knepley 
1339c1cad2e7SMatthew G. Knepley   ierr = DMPlexGetCone(dm, e, &cone);CHKERRQ(ierr);
1340c1cad2e7SMatthew G. Knepley   ierr = DMLabelGetValue(bodyLabel, cone[0], &bodyID_0);CHKERRQ(ierr);    // Do I need to check the other end?
1341c1cad2e7SMatthew G. Knepley   ierr = DMLabelGetValue(vertexLabel, cone[0], &vertexID_0);CHKERRQ(ierr);
1342c1cad2e7SMatthew G. Knepley   ierr = DMLabelGetValue(vertexLabel, cone[1], &vertexID_1);CHKERRQ(ierr);
1343c1cad2e7SMatthew G. Knepley   ierr = DMLabelGetValue(edgeLabel, cone[0], &edgeID_0);CHKERRQ(ierr);
1344c1cad2e7SMatthew G. Knepley   ierr = DMLabelGetValue(edgeLabel, cone[1], &edgeID_1);CHKERRQ(ierr);
1345c1cad2e7SMatthew G. Knepley   ierr = DMLabelGetValue(faceLabel, cone[0], &faceID_0);CHKERRQ(ierr);
1346c1cad2e7SMatthew G. Knepley   ierr = DMLabelGetValue(faceLabel, cone[1], &faceID_1);CHKERRQ(ierr);
1347c1cad2e7SMatthew G. Knepley 
1348c1cad2e7SMatthew G. Knepley   ierr = DMLabelSetValue(bodyLabel, e, bodyID_0);CHKERRQ(ierr);
1349c1cad2e7SMatthew G. Knepley 
1350c1cad2e7SMatthew G. Knepley   if (edgeID_0 == edgeID_1) { ierr = DMLabelSetValue(edgeLabel, e, edgeID_0);CHKERRQ(ierr); }
1351c1cad2e7SMatthew G. Knepley   else if (vertexID_0 > 0 && edgeID_1 > 0) { ierr = DMLabelSetValue(edgeLabel, e, edgeID_1);CHKERRQ(ierr); }
1352c1cad2e7SMatthew G. Knepley   else if (vertexID_1 > 0 && edgeID_0 > 0) { ierr = DMLabelSetValue(edgeLabel, e, edgeID_0);CHKERRQ(ierr); }
1353c1cad2e7SMatthew G. Knepley   else { /* Do Nothing */ }
1354c1cad2e7SMatthew G. Knepley   }
1355c1cad2e7SMatthew G. Knepley 
1356c1cad2e7SMatthew G. Knepley   /* Set DMLabels for Cells */
1357c1cad2e7SMatthew G. Knepley   for (int f = fStart; f < fEnd; ++f){
1358c1cad2e7SMatthew G. Knepley   int             edgeID_0;
1359c1cad2e7SMatthew G. Knepley   PetscInt        triBodyVal, triFaceVal;
1360c1cad2e7SMatthew G. Knepley   PetscHashIter   TFLiter, TBLiter;
1361c1cad2e7SMatthew G. Knepley   PetscBool       TFLfound, TBLfound;
1362c1cad2e7SMatthew G. Knepley 
1363c1cad2e7SMatthew G. Knepley     // Convert to Hash Table
1364c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIFind(triFaceIDLabelMap, f - fStart, &TFLiter, &TFLfound);CHKERRQ(ierr);
1365c1cad2e7SMatthew G. Knepley   if (!TFLfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in triFaceIDLabelMap", f);
1366c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIGet(triFaceIDLabelMap, f - fStart, &triFaceVal);CHKERRQ(ierr);
1367c1cad2e7SMatthew G. Knepley 
1368c1cad2e7SMatthew G. Knepley   ierr = PetscHMapIFind(triBodyIDLabelMap, f - fStart, &TBLiter, &TBLfound);CHKERRQ(ierr);
1369c1cad2e7SMatthew G. Knepley   if (!TBLfound) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in triBodyIDLabelMap", f);
1370c1cad2e7SMatthew G. Knepley     ierr = PetscHMapIGet(triBodyIDLabelMap, f - fStart, &triBodyVal);CHKERRQ(ierr);
1371c1cad2e7SMatthew G. Knepley 
1372c1cad2e7SMatthew G. Knepley   ierr = DMLabelSetValue(bodyLabel, f, triBodyVal);CHKERRQ(ierr);
1373c1cad2e7SMatthew G. Knepley   ierr = DMLabelSetValue(faceLabel, f, triFaceVal);CHKERRQ(ierr);
1374c1cad2e7SMatthew G. Knepley 
1375c1cad2e7SMatthew G. Knepley   /* Finish Labeling previously unlabeled DMPlex Edges - Assumes Triangular Cell (3 Edges Max) */
1376c1cad2e7SMatthew G. Knepley   ierr = DMPlexGetCone(dm, f, &cone);CHKERRQ(ierr);
1377c1cad2e7SMatthew G. Knepley 
1378c1cad2e7SMatthew G. Knepley   for (int jj = 0; jj < 3; ++jj) {
1379c1cad2e7SMatthew G. Knepley     ierr = DMLabelGetValue(edgeLabel, cone[jj], &edgeID_0);CHKERRQ(ierr);
1380c1cad2e7SMatthew G. Knepley 
1381c1cad2e7SMatthew G. Knepley     if (edgeID_0 < 0) {
1382c1cad2e7SMatthew G. Knepley     ierr = DMLabelSetValue(bodyLabel, cone[jj], triBodyVal);CHKERRQ(ierr);
1383c1cad2e7SMatthew G. Knepley       ierr = DMLabelSetValue(faceLabel, cone[jj], triFaceVal);CHKERRQ(ierr);
1384c1cad2e7SMatthew G. Knepley     }
1385c1cad2e7SMatthew G. Knepley   }
1386c1cad2e7SMatthew G. Knepley   }
1387c1cad2e7SMatthew G. Knepley 
1388c1cad2e7SMatthew G. Knepley   *newdm = dm;
1389c1cad2e7SMatthew G. Knepley   PetscFunctionReturn(0);
1390c1cad2e7SMatthew G. Knepley }
13917bee2925SMatthew Knepley #endif
13927bee2925SMatthew Knepley 
1393c1cad2e7SMatthew G. Knepley /*@
1394c1cad2e7SMatthew G. Knepley   DMPlexInflateToGeomModel - Snaps the vertex coordinates of a DMPlex object representing the mesh to its geometry if some vertices depart from the model. This usually happens with non-conforming refinement.
1395c1cad2e7SMatthew G. Knepley 
1396c1cad2e7SMatthew G. Knepley   Collective on dm
1397c1cad2e7SMatthew G. Knepley 
1398c1cad2e7SMatthew G. Knepley   Input Parameter:
1399c1cad2e7SMatthew G. Knepley . dm - The uninflated DM object representing the mesh
1400c1cad2e7SMatthew G. Knepley 
1401c1cad2e7SMatthew G. Knepley   Output Parameter:
1402c1cad2e7SMatthew G. Knepley . dm - The inflated DM object representing the mesh
1403c1cad2e7SMatthew G. Knepley 
1404c1cad2e7SMatthew G. Knepley   Level: intermediate
1405c1cad2e7SMatthew G. Knepley 
1406c1cad2e7SMatthew G. Knepley .seealso: DMPLEX, DMCreate(), DMPlexCreateEGADS()
1407c1cad2e7SMatthew G. Knepley @*/
1408c1cad2e7SMatthew G. Knepley PetscErrorCode DMPlexInflateToGeomModel(DM dm)
1409c1cad2e7SMatthew G. Knepley {
1410c1cad2e7SMatthew G. Knepley #if defined(PETSC_HAVE_EGADS)
1411c1cad2e7SMatthew G. Knepley   /* EGADS Variables */
1412c1cad2e7SMatthew G. Knepley   ego            model, geom, body, face, edge;
1413c1cad2e7SMatthew G. Knepley   ego           *bodies;
1414c1cad2e7SMatthew G. Knepley   int            Nb, oclass, mtype, *senses;
1415c1cad2e7SMatthew G. Knepley   double         result[3];
1416c1cad2e7SMatthew G. Knepley   /* PETSc Variables */
1417c1cad2e7SMatthew G. Knepley   DM             cdm;
1418c1cad2e7SMatthew G. Knepley   PetscContainer modelObj;
1419c1cad2e7SMatthew G. Knepley   DMLabel        bodyLabel, faceLabel, edgeLabel, vertexLabel;
1420c1cad2e7SMatthew G. Knepley   Vec            coordinates;
1421c1cad2e7SMatthew G. Knepley   PetscScalar   *coords;
1422c1cad2e7SMatthew G. Knepley   PetscInt       bodyID, faceID, edgeID, vertexID;
1423c1cad2e7SMatthew G. Knepley   PetscInt       cdim, d, vStart, vEnd, v;
1424c1cad2e7SMatthew G. Knepley   PetscErrorCode ierr;
1425c1cad2e7SMatthew G. Knepley #endif
1426c1cad2e7SMatthew G. Knepley 
1427c1cad2e7SMatthew G. Knepley   PetscFunctionBegin;
1428c1cad2e7SMatthew G. Knepley #if defined(PETSC_HAVE_EGADS)
1429c1cad2e7SMatthew G. Knepley   ierr = PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj);CHKERRQ(ierr);
1430c1cad2e7SMatthew G. Knepley   if (!modelObj) PetscFunctionReturn(0);
1431c1cad2e7SMatthew G. Knepley   ierr = DMGetCoordinateDim(dm, &cdim);CHKERRQ(ierr);
1432c1cad2e7SMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
1433c1cad2e7SMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
1434c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr);
1435c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Face ID", &faceLabel);CHKERRQ(ierr);
1436c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Edge ID", &edgeLabel);CHKERRQ(ierr);
1437c1cad2e7SMatthew G. Knepley   ierr = DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel);CHKERRQ(ierr);
1438c1cad2e7SMatthew G. Knepley 
1439c1cad2e7SMatthew G. Knepley   ierr = PetscContainerGetPointer(modelObj, (void **) &model);CHKERRQ(ierr);
1440c1cad2e7SMatthew G. Knepley   ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr);
1441c1cad2e7SMatthew G. Knepley 
1442c1cad2e7SMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1443c1cad2e7SMatthew G. Knepley   ierr = VecGetArrayWrite(coordinates, &coords);CHKERRQ(ierr);
1444c1cad2e7SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
1445c1cad2e7SMatthew G. Knepley     PetscScalar *vcoords;
1446c1cad2e7SMatthew G. Knepley 
1447c1cad2e7SMatthew G. Knepley     ierr = DMLabelGetValue(bodyLabel, v, &bodyID);CHKERRQ(ierr);
1448c1cad2e7SMatthew G. Knepley     ierr = DMLabelGetValue(faceLabel, v, &faceID);CHKERRQ(ierr);
1449c1cad2e7SMatthew G. Knepley     ierr = DMLabelGetValue(edgeLabel, v, &edgeID);CHKERRQ(ierr);
1450c1cad2e7SMatthew G. Knepley     ierr = DMLabelGetValue(vertexLabel, v, &vertexID);CHKERRQ(ierr);
1451c1cad2e7SMatthew G. Knepley 
1452c1cad2e7SMatthew G. Knepley     if (bodyID >= Nb) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %D is not in [0, %d)", bodyID, Nb);
1453c1cad2e7SMatthew G. Knepley     body = bodies[bodyID];
1454c1cad2e7SMatthew G. Knepley 
1455c1cad2e7SMatthew G. Knepley     ierr = DMPlexPointLocalRef(cdm, v, coords, (void *) &vcoords);CHKERRQ(ierr);
1456c1cad2e7SMatthew G. Knepley     if (edgeID > 0) {
1457c1cad2e7SMatthew G. Knepley       /* Snap to EDGE at nearest location */
1458c1cad2e7SMatthew G. Knepley       double params[1];
1459c1cad2e7SMatthew G. Knepley       ierr = EG_objectBodyTopo(body, EDGE, edgeID, &edge);CHKERRQ(ierr);
1460c1cad2e7SMatthew G. Knepley       ierr = EG_invEvaluate(edge, vcoords, params, result);CHKERRQ(ierr); // Get (x,y,z) of nearest point on EDGE
1461c1cad2e7SMatthew G. Knepley       for (d = 0; d < cdim; ++d) vcoords[d] = result[d];
1462c1cad2e7SMatthew G. Knepley     } else if (faceID > 0) {
1463c1cad2e7SMatthew G. Knepley       /* Snap to FACE at nearest location */
1464c1cad2e7SMatthew G. Knepley       double params[2];
1465c1cad2e7SMatthew G. Knepley       ierr = EG_objectBodyTopo(body, FACE, faceID, &face);CHKERRQ(ierr);
1466c1cad2e7SMatthew G. Knepley       ierr = EG_invEvaluate(face, vcoords, params, result);CHKERRQ(ierr); // Get (x,y,z) of nearest point on FACE
1467c1cad2e7SMatthew G. Knepley       for (d = 0; d < cdim; ++d) vcoords[d] = result[d];
1468c1cad2e7SMatthew G. Knepley     }
1469c1cad2e7SMatthew G. Knepley   }
1470c1cad2e7SMatthew G. Knepley   ierr = VecRestoreArrayWrite(coordinates, &coords);CHKERRQ(ierr);
1471c1cad2e7SMatthew G. Knepley   /* Clear out global coordinates */
1472c1cad2e7SMatthew G. Knepley   ierr = VecDestroy(&dm->coordinates);CHKERRQ(ierr);
1473c1cad2e7SMatthew G. Knepley #endif
1474c1cad2e7SMatthew G. Knepley   PetscFunctionReturn(0);
1475c1cad2e7SMatthew G. Knepley }
1476c1cad2e7SMatthew G. Knepley 
14777bee2925SMatthew Knepley /*@C
1478c1cad2e7SMatthew G. Knepley   DMPlexCreateEGADSFromFile - Create a DMPlex mesh from an EGADS, IGES, or STEP file.
14797bee2925SMatthew Knepley 
14807bee2925SMatthew Knepley   Collective
14817bee2925SMatthew Knepley 
14827bee2925SMatthew Knepley   Input Parameters:
14837bee2925SMatthew Knepley + comm     - The MPI communicator
1484c1cad2e7SMatthew G. Knepley - filename - The name of the EGADS, IGES, or STEP file
14857bee2925SMatthew Knepley 
14867bee2925SMatthew Knepley   Output Parameter:
14877bee2925SMatthew Knepley . dm       - The DM object representing the mesh
14887bee2925SMatthew Knepley 
14897bee2925SMatthew Knepley   Level: beginner
14907bee2925SMatthew Knepley 
1491c1cad2e7SMatthew G. Knepley .seealso: DMPLEX, DMCreate(), DMPlexCreateEGADS(), DMPlexCreateEGADSLiteFromFile()
14927bee2925SMatthew Knepley @*/
14937bee2925SMatthew Knepley PetscErrorCode DMPlexCreateEGADSFromFile(MPI_Comm comm, const char filename[], DM *dm)
14947bee2925SMatthew Knepley {
14957bee2925SMatthew Knepley   PetscMPIInt    rank;
14967bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS)
14977bee2925SMatthew Knepley   ego            context= NULL, model = NULL;
14987bee2925SMatthew Knepley #endif
1499c1cad2e7SMatthew G. Knepley   PetscBool      printModel = PETSC_FALSE, tessModel = PETSC_FALSE, newModel = PETSC_FALSE;
15007bee2925SMatthew Knepley   PetscErrorCode ierr;
15017bee2925SMatthew Knepley 
15027bee2925SMatthew Knepley   PetscFunctionBegin;
15037bee2925SMatthew Knepley   PetscValidCharPointer(filename, 2);
15047bee2925SMatthew Knepley   ierr = PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_print_model", &printModel, NULL);CHKERRQ(ierr);
1505c1cad2e7SMatthew G. Knepley   ierr = PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_tess_model", &tessModel, NULL);CHKERRQ(ierr);
1506c1cad2e7SMatthew G. Knepley   ierr = PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_new_model", &newModel, NULL);CHKERRQ(ierr);
150755b25c41SPierre Jolivet   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
15087bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS)
1509dd400576SPatrick Sanan   if (rank == 0) {
15107bee2925SMatthew Knepley 
15117bee2925SMatthew Knepley     ierr = EG_open(&context);CHKERRQ(ierr);
15127bee2925SMatthew Knepley     ierr = EG_loadModel(context, 0, filename, &model);CHKERRQ(ierr);
1513c1cad2e7SMatthew G. Knepley     if (printModel) {ierr = DMPlexEGADSPrintModel_Internal(model);CHKERRQ(ierr);}
15147bee2925SMatthew Knepley 
15157bee2925SMatthew Knepley   }
1516c1cad2e7SMatthew G. Knepley   if (tessModel)     {ierr = DMPlexCreateEGADS_Tess_Internal(comm, context, model, dm);CHKERRQ(ierr);}
1517c1cad2e7SMatthew G. Knepley   else if (newModel) {ierr = DMPlexCreateEGADS(comm, context, model, dm);CHKERRQ(ierr);}
1518c1cad2e7SMatthew G. Knepley   else               {ierr = DMPlexCreateEGADS_Internal(comm, context, model, dm);CHKERRQ(ierr);}
1519c03d1177SJose E. Roman   PetscFunctionReturn(0);
15207bee2925SMatthew Knepley #else
1521c1cad2e7SMatthew G. Knepley   SETERRQ(comm, PETSC_ERR_SUP, "This method requires EGADS support. Reconfigure using --download-egads");
15227bee2925SMatthew Knepley #endif
15237bee2925SMatthew Knepley }
1524