xref: /petsc/src/dm/impls/plex/plexegads.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
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 
25*9371c9d4SSatish Balay PetscErrorCode DMPlexSnapToGeomModel_EGADS_Internal(DM dm, PetscInt p, ego model, PetscInt bodyID, PetscInt faceID, PetscInt edgeID, const PetscScalar mcoords[], PetscScalar gcoords[]) {
26c1cad2e7SMatthew G. Knepley   DM           cdm;
27c1cad2e7SMatthew G. Knepley   ego         *bodies;
28c1cad2e7SMatthew G. Knepley   ego          geom, body, obj;
29a5b23f4aSJose E. Roman   /* result has to hold derivatives, along with the value */
30c1cad2e7SMatthew G. Knepley   double       params[3], result[18], paramsV[16 * 3], resultV[16 * 3], range[4];
31c1cad2e7SMatthew G. Knepley   int          Nb, oclass, mtype, *senses, peri;
32c1cad2e7SMatthew G. Knepley   Vec          coordinatesLocal;
33c1cad2e7SMatthew G. Knepley   PetscScalar *coords = NULL;
34c1cad2e7SMatthew G. Knepley   PetscInt     Nv, v, Np = 0, pm;
35c1cad2e7SMatthew G. Knepley   PetscInt     dE, d;
36c1cad2e7SMatthew G. Knepley 
37c1cad2e7SMatthew G. Knepley   PetscFunctionBeginHot;
389566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
399566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &dE));
409566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinatesLocal));
419566063dSJacob Faibussowitsch   PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
4263a3b9bcSJacob Faibussowitsch   PetscCheck(bodyID < Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %" PetscInt_FMT " is not in [0, %d)", bodyID, Nb);
43c1cad2e7SMatthew G. Knepley   body = bodies[bodyID];
44c1cad2e7SMatthew G. Knepley 
45*9371c9d4SSatish Balay   if (edgeID >= 0) {
46*9371c9d4SSatish Balay     PetscCall(EG_objectBodyTopo(body, EDGE, edgeID, &obj));
47*9371c9d4SSatish Balay     Np = 1;
48*9371c9d4SSatish Balay   } else if (faceID >= 0) {
49*9371c9d4SSatish Balay     PetscCall(EG_objectBodyTopo(body, FACE, faceID, &obj));
50*9371c9d4SSatish Balay     Np = 2;
51*9371c9d4SSatish Balay   } else {
52c1cad2e7SMatthew G. Knepley     for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
53c1cad2e7SMatthew G. Knepley     PetscFunctionReturn(0);
54c1cad2e7SMatthew G. Knepley   }
55c1cad2e7SMatthew G. Knepley 
56c1cad2e7SMatthew G. Knepley   /* Calculate parameters (t or u,v) for vertices */
579566063dSJacob Faibussowitsch   PetscCall(DMPlexVecGetClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords));
58c1cad2e7SMatthew G. Knepley   Nv /= dE;
59c1cad2e7SMatthew G. Knepley   if (Nv == 1) {
609566063dSJacob Faibussowitsch     PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords));
61c1cad2e7SMatthew G. Knepley     for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
62c1cad2e7SMatthew G. Knepley     PetscFunctionReturn(0);
63c1cad2e7SMatthew G. Knepley   }
6463a3b9bcSJacob Faibussowitsch   PetscCheck(Nv <= 16, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " coordinates associated to point %" PetscInt_FMT, Nv, p);
65c1cad2e7SMatthew G. Knepley 
66c1cad2e7SMatthew G. Knepley   /* Correct EGADSlite 2pi bug when calculating nearest point on Periodic Surfaces */
679566063dSJacob Faibussowitsch   PetscCall(EG_getRange(obj, range, &peri));
68c1cad2e7SMatthew G. Knepley   for (v = 0; v < Nv; ++v) {
699566063dSJacob Faibussowitsch     PetscCall(EG_invEvaluate(obj, &coords[v * dE], &paramsV[v * 3], &resultV[v * 3]));
70c1cad2e7SMatthew G. Knepley #if 1
71c1cad2e7SMatthew G. Knepley     if (peri > 0) {
72*9371c9d4SSatish Balay       if (paramsV[v * 3 + 0] + 1.e-4 < range[0]) {
73*9371c9d4SSatish Balay         paramsV[v * 3 + 0] += 2. * PETSC_PI;
74*9371c9d4SSatish Balay       } else if (paramsV[v * 3 + 0] - 1.e-4 > range[1]) {
75*9371c9d4SSatish Balay         paramsV[v * 3 + 0] -= 2. * PETSC_PI;
76*9371c9d4SSatish Balay       }
77c1cad2e7SMatthew G. Knepley     }
78c1cad2e7SMatthew G. Knepley     if (peri > 1) {
79*9371c9d4SSatish Balay       if (paramsV[v * 3 + 1] + 1.e-4 < range[2]) {
80*9371c9d4SSatish Balay         paramsV[v * 3 + 1] += 2. * PETSC_PI;
81*9371c9d4SSatish Balay       } else if (paramsV[v * 3 + 1] - 1.e-4 > range[3]) {
82*9371c9d4SSatish Balay         paramsV[v * 3 + 1] -= 2. * PETSC_PI;
83*9371c9d4SSatish Balay       }
84c1cad2e7SMatthew G. Knepley     }
85c1cad2e7SMatthew G. Knepley #endif
86c1cad2e7SMatthew G. Knepley   }
879566063dSJacob Faibussowitsch   PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords));
88c1cad2e7SMatthew G. Knepley   /* Calculate parameters (t or u,v) for new vertex at edge midpoint */
89c1cad2e7SMatthew G. Knepley   for (pm = 0; pm < Np; ++pm) {
90c1cad2e7SMatthew G. Knepley     params[pm] = 0.;
91c1cad2e7SMatthew G. Knepley     for (v = 0; v < Nv; ++v) params[pm] += paramsV[v * 3 + pm];
92c1cad2e7SMatthew G. Knepley     params[pm] /= Nv;
93c1cad2e7SMatthew G. Knepley   }
9463a3b9bcSJacob Faibussowitsch   PetscCheck(!(params[0] < range[0]) && !(params[0] > range[1]), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " had bad interpolation", p);
951dca8a05SBarry Smith   PetscCheck(Np <= 1 || (params[1] >= range[2] && params[1] <= range[3]), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " had bad interpolation", p);
96c1cad2e7SMatthew G. Knepley   /* Put coordinates for new vertex in result[] */
979566063dSJacob Faibussowitsch   PetscCall(EG_evaluate(obj, params, result));
98c1cad2e7SMatthew G. Knepley   for (d = 0; d < dE; ++d) gcoords[d] = result[d];
99c1cad2e7SMatthew G. Knepley   PetscFunctionReturn(0);
100c1cad2e7SMatthew G. Knepley }
101c1cad2e7SMatthew G. Knepley #endif
102c1cad2e7SMatthew G. Knepley 
103a8ededdfSMatthew G. Knepley /*@
104a8ededdfSMatthew 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.
105a8ededdfSMatthew G. Knepley 
106a8ededdfSMatthew G. Knepley   Not collective
107a8ededdfSMatthew G. Knepley 
108a8ededdfSMatthew G. Knepley   Input Parameters:
109a8ededdfSMatthew G. Knepley + dm      - The DMPlex object
110a8ededdfSMatthew G. Knepley . p       - The mesh point
111d410b0cfSMatthew G. Knepley . dE      - The coordinate dimension
112a8ededdfSMatthew G. Knepley - mcoords - A coordinate point lying on the mesh point
113a8ededdfSMatthew G. Knepley 
114a8ededdfSMatthew G. Knepley   Output Parameter:
115a8ededdfSMatthew G. Knepley . gcoords - The closest coordinate point on the geometry model associated with 'p' to the given point
116a8ededdfSMatthew G. Knepley 
117d410b0cfSMatthew 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.
118a8ededdfSMatthew G. Knepley 
119a8ededdfSMatthew G. Knepley   Level: intermediate
120a8ededdfSMatthew G. Knepley 
121db781477SPatrick Sanan .seealso: `DMRefine()`, `DMPlexCreate()`, `DMPlexSetRefinementUniform()`
122a8ededdfSMatthew G. Knepley @*/
123*9371c9d4SSatish Balay PetscErrorCode DMPlexSnapToGeomModel(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[]) {
124d410b0cfSMatthew G. Knepley   PetscInt d;
125a8ededdfSMatthew G. Knepley 
126c1cad2e7SMatthew G. Knepley   PetscFunctionBeginHot;
127a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS
128c1cad2e7SMatthew G. Knepley   {
129c1cad2e7SMatthew G. Knepley     DM_Plex       *plex = (DM_Plex *)dm->data;
130c1cad2e7SMatthew G. Knepley     DMLabel        bodyLabel, faceLabel, edgeLabel;
131c1cad2e7SMatthew G. Knepley     PetscInt       bodyID, faceID, edgeID;
132c1cad2e7SMatthew G. Knepley     PetscContainer modelObj;
133c1cad2e7SMatthew G. Knepley     ego            model;
134c1cad2e7SMatthew G. Knepley     PetscBool      islite = PETSC_FALSE;
135c1cad2e7SMatthew G. Knepley 
1369566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
1379566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
1389566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
139c1cad2e7SMatthew G. Knepley     if (!bodyLabel || !faceLabel || !edgeLabel || plex->ignoreModel) {
140f0fcf0adSMatthew G. Knepley       for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
141a8ededdfSMatthew G. Knepley       PetscFunctionReturn(0);
142a8ededdfSMatthew G. Knepley     }
1439566063dSJacob Faibussowitsch     PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
144c1cad2e7SMatthew G. Knepley     if (!modelObj) {
1459566063dSJacob Faibussowitsch       PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSLite Model", (PetscObject *)&modelObj));
146c1cad2e7SMatthew G. Knepley       islite = PETSC_TRUE;
147c1cad2e7SMatthew G. Knepley     }
148c1cad2e7SMatthew G. Knepley     if (!modelObj) {
149c1cad2e7SMatthew G. Knepley       for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
150c1cad2e7SMatthew G. Knepley       PetscFunctionReturn(0);
151c1cad2e7SMatthew G. Knepley     }
1529566063dSJacob Faibussowitsch     PetscCall(PetscContainerGetPointer(modelObj, (void **)&model));
1539566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(bodyLabel, p, &bodyID));
1549566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(faceLabel, p, &faceID));
1559566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(edgeLabel, p, &edgeID));
156c1cad2e7SMatthew G. Knepley     /* Allows for "Connective" Plex Edges present in models with multiple non-touching Entities */
157c1cad2e7SMatthew G. Knepley     if (bodyID < 0) {
158f0fcf0adSMatthew G. Knepley       for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
159f0fcf0adSMatthew G. Knepley       PetscFunctionReturn(0);
160a8ededdfSMatthew G. Knepley     }
1619566063dSJacob Faibussowitsch     if (islite) PetscCall(DMPlexSnapToGeomModel_EGADSLite_Internal(dm, p, model, bodyID, faceID, edgeID, mcoords, gcoords));
1629566063dSJacob Faibussowitsch     else PetscCall(DMPlexSnapToGeomModel_EGADS_Internal(dm, p, model, bodyID, faceID, edgeID, mcoords, gcoords));
163f0fcf0adSMatthew G. Knepley   }
164a8ededdfSMatthew G. Knepley #else
165f0fcf0adSMatthew G. Knepley   for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d];
166a8ededdfSMatthew G. Knepley #endif
167a8ededdfSMatthew G. Knepley   PetscFunctionReturn(0);
168a8ededdfSMatthew G. Knepley }
1697bee2925SMatthew Knepley 
1707bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS)
171*9371c9d4SSatish Balay static PetscErrorCode DMPlexEGADSPrintModel_Internal(ego model) {
1727bee2925SMatthew Knepley   ego geom, *bodies, *objs, *nobjs, *mobjs, *lobjs;
1737bee2925SMatthew Knepley   int oclass, mtype, *senses;
1747bee2925SMatthew Knepley   int Nb, b;
1757bee2925SMatthew Knepley 
1767bee2925SMatthew Knepley   PetscFunctionBeginUser;
1777bee2925SMatthew Knepley   /* test bodyTopo functions */
1789566063dSJacob Faibussowitsch   PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
1799566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of BODIES (nbodies): %d \n", Nb));
1807bee2925SMatthew Knepley 
1817bee2925SMatthew Knepley   for (b = 0; b < Nb; ++b) {
1827bee2925SMatthew Knepley     ego body = bodies[b];
1837bee2925SMatthew Knepley     int id, Nsh, Nf, Nl, l, Ne, e, Nv, v;
1847bee2925SMatthew Knepley 
1857bee2925SMatthew Knepley     /* Output Basic Model Topology */
1869566063dSJacob Faibussowitsch     PetscCall(EG_getBodyTopos(body, NULL, SHELL, &Nsh, &objs));
1879566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "   Number of SHELLS: %d \n", Nsh));
1887bee2925SMatthew Knepley     EG_free(objs);
1897bee2925SMatthew Knepley 
1909566063dSJacob Faibussowitsch     PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &objs));
1919566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "   Number of FACES: %d \n", Nf));
1927bee2925SMatthew Knepley     EG_free(objs);
1937bee2925SMatthew Knepley 
1949566063dSJacob Faibussowitsch     PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
1959566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "   Number of LOOPS: %d \n", Nl));
1967bee2925SMatthew Knepley 
1979566063dSJacob Faibussowitsch     PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &objs));
1989566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "   Number of EDGES: %d \n", Ne));
1997bee2925SMatthew Knepley     EG_free(objs);
2007bee2925SMatthew Knepley 
2019566063dSJacob Faibussowitsch     PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &objs));
2029566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF, "   Number of NODES: %d \n", Nv));
2037bee2925SMatthew Knepley     EG_free(objs);
2047bee2925SMatthew Knepley 
2057bee2925SMatthew Knepley     for (l = 0; l < Nl; ++l) {
2067bee2925SMatthew Knepley       ego loop = lobjs[l];
2077bee2925SMatthew Knepley 
2087bee2925SMatthew Knepley       id = EG_indexBodyTopo(body, loop);
2099566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF, "          LOOP ID: %d\n", id));
2107bee2925SMatthew Knepley 
2117bee2925SMatthew Knepley       /* Get EDGE info which associated with the current LOOP */
2129566063dSJacob Faibussowitsch       PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
2137bee2925SMatthew Knepley 
2147bee2925SMatthew Knepley       for (e = 0; e < Ne; ++e) {
2157bee2925SMatthew Knepley         ego    edge      = objs[e];
2167bee2925SMatthew Knepley         double range[4]  = {0., 0., 0., 0.};
2177bee2925SMatthew Knepley         double point[3]  = {0., 0., 0.};
218266cfabeSMatthew G. Knepley         double params[3] = {0., 0., 0.};
219266cfabeSMatthew G. Knepley         double result[18];
2207bee2925SMatthew Knepley         int    peri;
2217bee2925SMatthew Knepley 
2225f80ce2aSJacob Faibussowitsch         id = EG_indexBodyTopo(body, edge);
2239566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "            EDGE ID: %d (%d)\n", id, e));
2247bee2925SMatthew Knepley 
2259566063dSJacob Faibussowitsch         PetscCall(EG_getRange(edge, range, &peri));
2269566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_SELF, "  Range = %lf, %lf, %lf, %lf \n", range[0], range[1], range[2], range[3]));
2277bee2925SMatthew Knepley 
2287bee2925SMatthew Knepley         /* Get NODE info which associated with the current EDGE */
2299566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
230266cfabeSMatthew G. Knepley         if (mtype == DEGENERATE) {
2319566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "  EDGE %d is DEGENERATE \n", id));
232266cfabeSMatthew G. Knepley         } else {
233266cfabeSMatthew G. Knepley           params[0] = range[0];
2349566063dSJacob Faibussowitsch           PetscCall(EG_evaluate(edge, params, result));
2359566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "   between (%lf, %lf, %lf)", result[0], result[1], result[2]));
236266cfabeSMatthew G. Knepley           params[0] = range[1];
2379566063dSJacob Faibussowitsch           PetscCall(EG_evaluate(edge, params, result));
2389566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, " and (%lf, %lf, %lf)\n", result[0], result[1], result[2]));
239266cfabeSMatthew G. Knepley         }
2407bee2925SMatthew Knepley 
2417bee2925SMatthew Knepley         for (v = 0; v < Nv; ++v) {
2427bee2925SMatthew Knepley           ego    vertex = nobjs[v];
2437bee2925SMatthew Knepley           double limits[4];
2447bee2925SMatthew Knepley           int    dummy;
2457bee2925SMatthew Knepley 
2469566063dSJacob Faibussowitsch           PetscCall(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses));
2477bee2925SMatthew Knepley           id = EG_indexBodyTopo(body, vertex);
2489566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "              NODE ID: %d \n", id));
2499566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF, "                 (x, y, z) = (%lf, %lf, %lf) \n", limits[0], limits[1], limits[2]));
2507bee2925SMatthew Knepley 
2517bee2925SMatthew Knepley           point[0] = point[0] + limits[0];
2527bee2925SMatthew Knepley           point[1] = point[1] + limits[1];
2537bee2925SMatthew Knepley           point[2] = point[2] + limits[2];
2547bee2925SMatthew Knepley         }
2557bee2925SMatthew Knepley       }
2567bee2925SMatthew Knepley     }
2577bee2925SMatthew Knepley     EG_free(lobjs);
2587bee2925SMatthew Knepley   }
2597bee2925SMatthew Knepley   PetscFunctionReturn(0);
2607bee2925SMatthew Knepley }
2617bee2925SMatthew Knepley 
262*9371c9d4SSatish Balay static PetscErrorCode DMPlexEGADSDestroy_Private(void *context) {
2637bee2925SMatthew Knepley   if (context) EG_close((ego)context);
2647bee2925SMatthew Knepley   return 0;
2657bee2925SMatthew Knepley }
2667bee2925SMatthew Knepley 
267*9371c9d4SSatish Balay static PetscErrorCode DMPlexCreateEGADS_Internal(MPI_Comm comm, ego context, ego model, DM *newdm) {
268c1cad2e7SMatthew G. Knepley   DMLabel     bodyLabel, faceLabel, edgeLabel, vertexLabel;
2697bee2925SMatthew Knepley   PetscInt    cStart, cEnd, c;
2707bee2925SMatthew Knepley   /* EGADSLite variables */
2717bee2925SMatthew Knepley   ego         geom, *bodies, *objs, *nobjs, *mobjs, *lobjs;
2727bee2925SMatthew Knepley   int         oclass, mtype, nbodies, *senses;
2737bee2925SMatthew Knepley   int         b;
2747bee2925SMatthew Knepley   /* PETSc variables */
2757bee2925SMatthew Knepley   DM          dm;
2767bee2925SMatthew Knepley   PetscHMapI  edgeMap = NULL;
2777bee2925SMatthew 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;
2787bee2925SMatthew Knepley   PetscInt   *cells = NULL, *cone = NULL;
2797bee2925SMatthew Knepley   PetscReal  *coords = NULL;
2807bee2925SMatthew Knepley   PetscMPIInt rank;
2817bee2925SMatthew Knepley 
2827bee2925SMatthew Knepley   PetscFunctionBegin;
2839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
284dd400576SPatrick Sanan   if (rank == 0) {
285266cfabeSMatthew G. Knepley     const PetscInt debug = 0;
286266cfabeSMatthew G. Knepley 
2877bee2925SMatthew Knepley     /* ---------------------------------------------------------------------------------------------------
2887bee2925SMatthew Knepley     Generate Petsc Plex
2897bee2925SMatthew Knepley       Get all Nodes in model, record coordinates in a correctly formatted array
2907bee2925SMatthew Knepley       Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array
2917bee2925SMatthew Knepley       We need to uniformly refine the initial geometry to guarantee a valid mesh
2927bee2925SMatthew Knepley     */
2937bee2925SMatthew Knepley 
2947bee2925SMatthew Knepley     /* Calculate cell and vertex sizes */
2959566063dSJacob Faibussowitsch     PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses));
2969566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&edgeMap));
2977bee2925SMatthew Knepley     numEdges = 0;
2987bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
2997bee2925SMatthew Knepley       ego body = bodies[b];
3007bee2925SMatthew Knepley       int id, Nl, l, Nv, v;
3017bee2925SMatthew Knepley 
3029566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
3037bee2925SMatthew Knepley       for (l = 0; l < Nl; ++l) {
3047bee2925SMatthew Knepley         ego loop = lobjs[l];
3057bee2925SMatthew Knepley         int Ner  = 0, Ne, e, Nc;
3067bee2925SMatthew Knepley 
3079566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
3087bee2925SMatthew Knepley         for (e = 0; e < Ne; ++e) {
3097bee2925SMatthew Knepley           ego           edge = objs[e];
3107bee2925SMatthew Knepley           int           Nv, id;
3117bee2925SMatthew Knepley           PetscHashIter iter;
3127bee2925SMatthew Knepley           PetscBool     found;
3137bee2925SMatthew Knepley 
3149566063dSJacob Faibussowitsch           PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
3157bee2925SMatthew Knepley           if (mtype == DEGENERATE) continue;
3165f80ce2aSJacob Faibussowitsch           id = EG_indexBodyTopo(body, edge);
3179566063dSJacob Faibussowitsch           PetscCall(PetscHMapIFind(edgeMap, id - 1, &iter, &found));
3189566063dSJacob Faibussowitsch           if (!found) PetscCall(PetscHMapISet(edgeMap, id - 1, numEdges++));
3197bee2925SMatthew Knepley           ++Ner;
3207bee2925SMatthew Knepley         }
321*9371c9d4SSatish Balay         if (Ner == 2) {
322*9371c9d4SSatish Balay           Nc = 2;
323*9371c9d4SSatish Balay         } else if (Ner == 3) {
324*9371c9d4SSatish Balay           Nc = 4;
325*9371c9d4SSatish Balay         } else if (Ner == 4) {
326*9371c9d4SSatish Balay           Nc = 8;
327*9371c9d4SSatish Balay           ++numQuads;
328*9371c9d4SSatish Balay         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot support loop with %d edges", Ner);
3297bee2925SMatthew Knepley         numCells += Nc;
3307bee2925SMatthew Knepley         newCells += Nc - 1;
3317bee2925SMatthew Knepley         maxCorners = PetscMax(Ner * 2 + 1, maxCorners);
3327bee2925SMatthew Knepley       }
3339566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
3347bee2925SMatthew Knepley       for (v = 0; v < Nv; ++v) {
3357bee2925SMatthew Knepley         ego vertex = nobjs[v];
3367bee2925SMatthew Knepley 
3377bee2925SMatthew Knepley         id          = EG_indexBodyTopo(body, vertex);
3387bee2925SMatthew Knepley         /* TODO: Instead of assuming contiguous ids, we could use a hash table */
3397bee2925SMatthew Knepley         numVertices = PetscMax(id, numVertices);
3407bee2925SMatthew Knepley       }
3417bee2925SMatthew Knepley       EG_free(lobjs);
3427bee2925SMatthew Knepley       EG_free(nobjs);
3437bee2925SMatthew Knepley     }
3449566063dSJacob Faibussowitsch     PetscCall(PetscHMapIGetSize(edgeMap, &numEdges));
3457bee2925SMatthew Knepley     newVertices = numEdges + numQuads;
3467bee2925SMatthew Knepley     numVertices += newVertices;
3477bee2925SMatthew Knepley 
3487bee2925SMatthew Knepley     dim        = 2; /* Assume 3D Models :: Need to update to handle 2D Models in the future */
3497bee2925SMatthew Knepley     cdim       = 3; /* Assume 3D Models :: Need to update to handle 2D Models in the future */
3507bee2925SMatthew Knepley     numCorners = 3; /* Split cells into triangles */
3519566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(numVertices * cdim, &coords, numCells * numCorners, &cells, maxCorners, &cone));
3527bee2925SMatthew Knepley 
3537bee2925SMatthew Knepley     /* Get vertex coordinates */
3547bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
3557bee2925SMatthew Knepley       ego body = bodies[b];
3567bee2925SMatthew Knepley       int id, Nv, v;
3577bee2925SMatthew Knepley 
3589566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
3597bee2925SMatthew Knepley       for (v = 0; v < Nv; ++v) {
3607bee2925SMatthew Knepley         ego    vertex = nobjs[v];
3617bee2925SMatthew Knepley         double limits[4];
3627bee2925SMatthew Knepley         int    dummy;
3637bee2925SMatthew Knepley 
3649566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses));
3655f80ce2aSJacob Faibussowitsch         id                          = EG_indexBodyTopo(body, vertex);
3667bee2925SMatthew Knepley         coords[(id - 1) * cdim + 0] = limits[0];
3677bee2925SMatthew Knepley         coords[(id - 1) * cdim + 1] = limits[1];
3687bee2925SMatthew Knepley         coords[(id - 1) * cdim + 2] = limits[2];
3697bee2925SMatthew Knepley       }
3707bee2925SMatthew Knepley       EG_free(nobjs);
3717bee2925SMatthew Knepley     }
3729566063dSJacob Faibussowitsch     PetscCall(PetscHMapIClear(edgeMap));
3737bee2925SMatthew Knepley     fOff     = numVertices - newVertices + numEdges;
3747bee2925SMatthew Knepley     numEdges = 0;
3757bee2925SMatthew Knepley     numQuads = 0;
3767bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
3777bee2925SMatthew Knepley       ego body = bodies[b];
3787bee2925SMatthew Knepley       int Nl, l;
3797bee2925SMatthew Knepley 
3809566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
3817bee2925SMatthew Knepley       for (l = 0; l < Nl; ++l) {
3827bee2925SMatthew Knepley         ego loop = lobjs[l];
3837bee2925SMatthew Knepley         int lid, Ner = 0, Ne, e;
3847bee2925SMatthew Knepley 
3855f80ce2aSJacob Faibussowitsch         lid = EG_indexBodyTopo(body, loop);
3869566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
3877bee2925SMatthew Knepley         for (e = 0; e < Ne; ++e) {
3887bee2925SMatthew Knepley           ego           edge = objs[e];
3897bee2925SMatthew Knepley           int           eid, Nv;
3907bee2925SMatthew Knepley           PetscHashIter iter;
3917bee2925SMatthew Knepley           PetscBool     found;
3927bee2925SMatthew Knepley 
3939566063dSJacob Faibussowitsch           PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
3947bee2925SMatthew Knepley           if (mtype == DEGENERATE) continue;
3957bee2925SMatthew Knepley           ++Ner;
3965f80ce2aSJacob Faibussowitsch           eid = EG_indexBodyTopo(body, edge);
3979566063dSJacob Faibussowitsch           PetscCall(PetscHMapIFind(edgeMap, eid - 1, &iter, &found));
3987bee2925SMatthew Knepley           if (!found) {
3997bee2925SMatthew Knepley             PetscInt v = numVertices - newVertices + numEdges;
4007bee2925SMatthew Knepley             double   range[4], params[3] = {0., 0., 0.}, result[18];
4017bee2925SMatthew Knepley             int      periodic[2];
4027bee2925SMatthew Knepley 
4039566063dSJacob Faibussowitsch             PetscCall(PetscHMapISet(edgeMap, eid - 1, numEdges++));
4049566063dSJacob Faibussowitsch             PetscCall(EG_getRange(edge, range, periodic));
4057bee2925SMatthew Knepley             params[0] = 0.5 * (range[0] + range[1]);
4069566063dSJacob Faibussowitsch             PetscCall(EG_evaluate(edge, params, result));
4077bee2925SMatthew Knepley             coords[v * cdim + 0] = result[0];
4087bee2925SMatthew Knepley             coords[v * cdim + 1] = result[1];
4097bee2925SMatthew Knepley             coords[v * cdim + 2] = result[2];
4107bee2925SMatthew Knepley           }
4117bee2925SMatthew Knepley         }
4127bee2925SMatthew Knepley         if (Ner == 4) {
4137bee2925SMatthew Knepley           PetscInt v = fOff + numQuads++;
414266cfabeSMatthew G. Knepley           ego     *fobjs, face;
4157bee2925SMatthew Knepley           double   range[4], params[3] = {0., 0., 0.}, result[18];
416266cfabeSMatthew G. Knepley           int      Nf, fid, periodic[2];
4177bee2925SMatthew Knepley 
4189566063dSJacob Faibussowitsch           PetscCall(EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs));
419266cfabeSMatthew G. Knepley           face = fobjs[0];
4205f80ce2aSJacob Faibussowitsch           fid  = EG_indexBodyTopo(body, face);
42108401ef6SPierre Jolivet           PetscCheck(Nf == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Loop %d has %d faces, instead of 1 (%d)", lid - 1, Nf, fid);
4229566063dSJacob Faibussowitsch           PetscCall(EG_getRange(face, range, periodic));
4237bee2925SMatthew Knepley           params[0] = 0.5 * (range[0] + range[1]);
4247bee2925SMatthew Knepley           params[1] = 0.5 * (range[2] + range[3]);
4259566063dSJacob Faibussowitsch           PetscCall(EG_evaluate(face, params, result));
4267bee2925SMatthew Knepley           coords[v * cdim + 0] = result[0];
4277bee2925SMatthew Knepley           coords[v * cdim + 1] = result[1];
4287bee2925SMatthew Knepley           coords[v * cdim + 2] = result[2];
4297bee2925SMatthew Knepley         }
4307bee2925SMatthew Knepley       }
4317bee2925SMatthew Knepley     }
4321dca8a05SBarry Smith     PetscCheck(numEdges + numQuads == newVertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of new vertices %" PetscInt_FMT " != %" PetscInt_FMT " previous count", numEdges + numQuads, newVertices);
4337bee2925SMatthew Knepley 
4347bee2925SMatthew Knepley     /* Get cell vertices by traversing loops */
4357bee2925SMatthew Knepley     numQuads = 0;
4367bee2925SMatthew Knepley     cOff     = 0;
4377bee2925SMatthew Knepley     for (b = 0; b < nbodies; ++b) {
4387bee2925SMatthew Knepley       ego body = bodies[b];
4397bee2925SMatthew Knepley       int id, Nl, l;
4407bee2925SMatthew Knepley 
4419566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
4427bee2925SMatthew Knepley       for (l = 0; l < Nl; ++l) {
4437bee2925SMatthew Knepley         ego loop = lobjs[l];
4447bee2925SMatthew Knepley         int lid, Ner = 0, Ne, e, nc = 0, c, Nt, t;
4457bee2925SMatthew Knepley 
4465f80ce2aSJacob Faibussowitsch         lid = EG_indexBodyTopo(body, loop);
4479566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
4487bee2925SMatthew Knepley 
4497bee2925SMatthew Knepley         for (e = 0; e < Ne; ++e) {
4507bee2925SMatthew Knepley           ego edge = objs[e];
4517bee2925SMatthew Knepley           int points[3];
4527bee2925SMatthew Knepley           int eid, Nv, v, tmp;
4537bee2925SMatthew Knepley 
4547bee2925SMatthew Knepley           eid = EG_indexBodyTopo(body, edge);
4559566063dSJacob Faibussowitsch           PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
456266cfabeSMatthew G. Knepley           if (mtype == DEGENERATE) continue;
457266cfabeSMatthew G. Knepley           else ++Ner;
45808401ef6SPierre Jolivet           PetscCheck(Nv == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Edge %d has %d vertices != 2", eid, Nv);
4597bee2925SMatthew Knepley 
4607bee2925SMatthew Knepley           for (v = 0; v < Nv; ++v) {
4617bee2925SMatthew Knepley             ego vertex = nobjs[v];
4627bee2925SMatthew Knepley 
4637bee2925SMatthew Knepley             id            = EG_indexBodyTopo(body, vertex);
4647bee2925SMatthew Knepley             points[v * 2] = id - 1;
4657bee2925SMatthew Knepley           }
4667bee2925SMatthew Knepley           {
4677bee2925SMatthew Knepley             PetscInt edgeNum;
4687bee2925SMatthew Knepley 
4699566063dSJacob Faibussowitsch             PetscCall(PetscHMapIGet(edgeMap, eid - 1, &edgeNum));
4707bee2925SMatthew Knepley             points[1] = numVertices - newVertices + edgeNum;
4717bee2925SMatthew Knepley           }
4727bee2925SMatthew Knepley           /* EGADS loops are not oriented, but seem to be in order, so we must piece them together */
4737bee2925SMatthew Knepley           if (!nc) {
4747bee2925SMatthew Knepley             for (v = 0; v < Nv + 1; ++v) cone[nc++] = points[v];
4757bee2925SMatthew Knepley           } else {
476*9371c9d4SSatish Balay             if (cone[nc - 1] == points[0]) {
477*9371c9d4SSatish Balay               cone[nc++] = points[1];
478*9371c9d4SSatish Balay               if (cone[0] != points[2]) cone[nc++] = points[2];
479*9371c9d4SSatish Balay             } else if (cone[nc - 1] == points[2]) {
480*9371c9d4SSatish Balay               cone[nc++] = points[1];
481*9371c9d4SSatish Balay               if (cone[0] != points[0]) cone[nc++] = points[0];
482*9371c9d4SSatish Balay             } else if (cone[nc - 3] == points[0]) {
483*9371c9d4SSatish Balay               tmp          = cone[nc - 3];
484*9371c9d4SSatish Balay               cone[nc - 3] = cone[nc - 1];
485*9371c9d4SSatish Balay               cone[nc - 1] = tmp;
486*9371c9d4SSatish Balay               cone[nc++]   = points[1];
487*9371c9d4SSatish Balay               if (cone[0] != points[2]) cone[nc++] = points[2];
488*9371c9d4SSatish Balay             } else if (cone[nc - 3] == points[2]) {
489*9371c9d4SSatish Balay               tmp          = cone[nc - 3];
490*9371c9d4SSatish Balay               cone[nc - 3] = cone[nc - 1];
491*9371c9d4SSatish Balay               cone[nc - 1] = tmp;
492*9371c9d4SSatish Balay               cone[nc++]   = points[1];
493*9371c9d4SSatish Balay               if (cone[0] != points[0]) cone[nc++] = points[0];
494*9371c9d4SSatish Balay             } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %d does not match its predecessor", eid);
4957bee2925SMatthew Knepley           }
4967bee2925SMatthew Knepley         }
49763a3b9bcSJacob Faibussowitsch         PetscCheck(nc == 2 * Ner, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %" PetscInt_FMT " != %" PetscInt_FMT, nc, 2 * Ner);
4987bee2925SMatthew Knepley         if (Ner == 4) { cone[nc++] = numVertices - newVertices + numEdges + numQuads++; }
49963a3b9bcSJacob Faibussowitsch         PetscCheck(nc <= maxCorners, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %" PetscInt_FMT " > %" PetscInt_FMT " max", nc, maxCorners);
5007bee2925SMatthew Knepley         /* Triangulate the loop */
5017bee2925SMatthew Knepley         switch (Ner) {
5027bee2925SMatthew Knepley         case 2: /* Bi-Segment -> 2 triangles */
5037bee2925SMatthew Knepley           Nt                           = 2;
5047bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[0];
5057bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[1];
5067bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[2];
5077bee2925SMatthew Knepley           ++cOff;
5087bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[0];
5097bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[2];
5107bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[3];
5117bee2925SMatthew Knepley           ++cOff;
5127bee2925SMatthew Knepley           break;
5137bee2925SMatthew Knepley         case 3: /* Triangle   -> 4 triangles */
5147bee2925SMatthew Knepley           Nt                           = 4;
5157bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[0];
5167bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[1];
5177bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[5];
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[5];
5247bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[3];
5257bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[4];
5267bee2925SMatthew Knepley           ++cOff;
5277bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[1];
5287bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[3];
5297bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[5];
5307bee2925SMatthew Knepley           ++cOff;
5317bee2925SMatthew Knepley           break;
5327bee2925SMatthew Knepley         case 4: /* Quad       -> 8 triangles */
5337bee2925SMatthew Knepley           Nt                           = 8;
5347bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[0];
5357bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[1];
5367bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[7];
5377bee2925SMatthew Knepley           ++cOff;
5387bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[1];
5397bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[2];
5407bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[3];
5417bee2925SMatthew Knepley           ++cOff;
5427bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[3];
5437bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[4];
5447bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[5];
5457bee2925SMatthew Knepley           ++cOff;
5467bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[5];
5477bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[6];
5487bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[7];
5497bee2925SMatthew Knepley           ++cOff;
5507bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[8];
5517bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[1];
5527bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[3];
5537bee2925SMatthew Knepley           ++cOff;
5547bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[8];
5557bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[3];
5567bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[5];
5577bee2925SMatthew Knepley           ++cOff;
5587bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[8];
5597bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[5];
5607bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[7];
5617bee2925SMatthew Knepley           ++cOff;
5627bee2925SMatthew Knepley           cells[cOff * numCorners + 0] = cone[8];
5637bee2925SMatthew Knepley           cells[cOff * numCorners + 1] = cone[7];
5647bee2925SMatthew Knepley           cells[cOff * numCorners + 2] = cone[1];
5657bee2925SMatthew Knepley           ++cOff;
5667bee2925SMatthew Knepley           break;
56798921bdaSJacob Faibussowitsch         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d edges, which we do not support", lid, Ner);
5687bee2925SMatthew Knepley         }
569266cfabeSMatthew G. Knepley         if (debug) {
5707bee2925SMatthew Knepley           for (t = 0; t < Nt; ++t) {
57163a3b9bcSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, "  LOOP Corner NODEs Triangle %" PetscInt_FMT " (", t));
5727bee2925SMatthew Knepley             for (c = 0; c < numCorners; ++c) {
5739566063dSJacob Faibussowitsch               if (c > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", "));
57463a3b9bcSJacob Faibussowitsch               PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT, cells[(cOff - Nt + t) * numCorners + c]));
5757bee2925SMatthew Knepley             }
5769566063dSJacob Faibussowitsch             PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n"));
5777bee2925SMatthew Knepley           }
5787bee2925SMatthew Knepley         }
579266cfabeSMatthew G. Knepley       }
5807bee2925SMatthew Knepley       EG_free(lobjs);
5817bee2925SMatthew Knepley     }
5827bee2925SMatthew Knepley   }
58363a3b9bcSJacob Faibussowitsch   PetscCheck(cOff == numCells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count of total cells %" PetscInt_FMT " != %" PetscInt_FMT " previous count", cOff, numCells);
5849566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numVertices, numCorners, PETSC_TRUE, cells, cdim, coords, &dm));
5859566063dSJacob Faibussowitsch   PetscCall(PetscFree3(coords, cells, cone));
58663a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(dm, " Total Number of Unique Cells    = %" PetscInt_FMT " (%" PetscInt_FMT ")\n", numCells, newCells));
58763a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(dm, " Total Number of Unique Vertices = %" PetscInt_FMT " (%" PetscInt_FMT ")\n", numVertices, newVertices));
5887bee2925SMatthew Knepley   /* Embed EGADS model in DM */
5897bee2925SMatthew Knepley   {
5907bee2925SMatthew Knepley     PetscContainer modelObj, contextObj;
5917bee2925SMatthew Knepley 
5929566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &modelObj));
5939566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(modelObj, model));
5949566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Model", (PetscObject)modelObj));
5959566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&modelObj));
5967bee2925SMatthew Knepley 
5979566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &contextObj));
5989566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(contextObj, context));
5999566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private));
6009566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Context", (PetscObject)contextObj));
6019566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&contextObj));
6027bee2925SMatthew Knepley   }
6037bee2925SMatthew Knepley   /* Label points */
6049566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Body ID"));
6059566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
6069566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Face ID"));
6079566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
6089566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Edge ID"));
6099566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
6109566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Vertex ID"));
6119566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
6127bee2925SMatthew Knepley   cOff = 0;
6137bee2925SMatthew Knepley   for (b = 0; b < nbodies; ++b) {
6147bee2925SMatthew Knepley     ego body = bodies[b];
6157bee2925SMatthew Knepley     int id, Nl, l;
6167bee2925SMatthew Knepley 
6179566063dSJacob Faibussowitsch     PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs));
6187bee2925SMatthew Knepley     for (l = 0; l < Nl; ++l) {
6197bee2925SMatthew Knepley       ego  loop = lobjs[l];
6207bee2925SMatthew Knepley       ego *fobjs;
6217bee2925SMatthew Knepley       int  lid, Nf, fid, Ner = 0, Ne, e, Nt = 0, t;
6227bee2925SMatthew Knepley 
6235f80ce2aSJacob Faibussowitsch       lid = EG_indexBodyTopo(body, loop);
6249566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs));
62508401ef6SPierre Jolivet       PetscCheck(Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf);
6265f80ce2aSJacob Faibussowitsch       fid = EG_indexBodyTopo(body, fobjs[0]);
6277bee2925SMatthew Knepley       EG_free(fobjs);
6289566063dSJacob Faibussowitsch       PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses));
6297bee2925SMatthew Knepley       for (e = 0; e < Ne; ++e) {
6307bee2925SMatthew Knepley         ego             edge = objs[e];
6317bee2925SMatthew Knepley         int             eid, Nv, v;
6327bee2925SMatthew Knepley         PetscInt        points[3], support[2], numEdges, edgeNum;
6337bee2925SMatthew Knepley         const PetscInt *edges;
6347bee2925SMatthew Knepley 
6357bee2925SMatthew Knepley         eid = EG_indexBodyTopo(body, edge);
6369566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
6377bee2925SMatthew Knepley         if (mtype == DEGENERATE) continue;
6387bee2925SMatthew Knepley         else ++Ner;
6397bee2925SMatthew Knepley         for (v = 0; v < Nv; ++v) {
6407bee2925SMatthew Knepley           ego vertex = nobjs[v];
6417bee2925SMatthew Knepley 
6427bee2925SMatthew Knepley           id = EG_indexBodyTopo(body, vertex);
6439566063dSJacob Faibussowitsch           PetscCall(DMLabelSetValue(edgeLabel, numCells + id - 1, eid));
6447bee2925SMatthew Knepley           points[v * 2] = numCells + id - 1;
6457bee2925SMatthew Knepley         }
6469566063dSJacob Faibussowitsch         PetscCall(PetscHMapIGet(edgeMap, eid - 1, &edgeNum));
6477bee2925SMatthew Knepley         points[1] = numCells + numVertices - newVertices + edgeNum;
6487bee2925SMatthew Knepley 
6499566063dSJacob Faibussowitsch         PetscCall(DMLabelSetValue(edgeLabel, points[1], eid));
6507bee2925SMatthew Knepley         support[0] = points[0];
6517bee2925SMatthew Knepley         support[1] = points[1];
6529566063dSJacob Faibussowitsch         PetscCall(DMPlexGetJoin(dm, 2, support, &numEdges, &edges));
65363a3b9bcSJacob Faibussowitsch         PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vertices (%" PetscInt_FMT ", %" PetscInt_FMT ") should only bound 1 edge, not %" PetscInt_FMT, support[0], support[1], numEdges);
6549566063dSJacob Faibussowitsch         PetscCall(DMLabelSetValue(edgeLabel, edges[0], eid));
6559566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges));
6567bee2925SMatthew Knepley         support[0] = points[1];
6577bee2925SMatthew Knepley         support[1] = points[2];
6589566063dSJacob Faibussowitsch         PetscCall(DMPlexGetJoin(dm, 2, support, &numEdges, &edges));
65963a3b9bcSJacob Faibussowitsch         PetscCheck(numEdges == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vertices (%" PetscInt_FMT ", %" PetscInt_FMT ") should only bound 1 edge, not %" PetscInt_FMT, support[0], support[1], numEdges);
6609566063dSJacob Faibussowitsch         PetscCall(DMLabelSetValue(edgeLabel, edges[0], eid));
6619566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges));
6627bee2925SMatthew Knepley       }
6637bee2925SMatthew Knepley       switch (Ner) {
6647bee2925SMatthew Knepley       case 2: Nt = 2; break;
6657bee2925SMatthew Knepley       case 3: Nt = 4; break;
6667bee2925SMatthew Knepley       case 4: Nt = 8; break;
66798921bdaSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Loop with %d edges is unsupported", Ner);
6687bee2925SMatthew Knepley       }
6697bee2925SMatthew Knepley       for (t = 0; t < Nt; ++t) {
6709566063dSJacob Faibussowitsch         PetscCall(DMLabelSetValue(bodyLabel, cOff + t, b));
6719566063dSJacob Faibussowitsch         PetscCall(DMLabelSetValue(faceLabel, cOff + t, fid));
6727bee2925SMatthew Knepley       }
6737bee2925SMatthew Knepley       cOff += Nt;
6747bee2925SMatthew Knepley     }
6757bee2925SMatthew Knepley     EG_free(lobjs);
6767bee2925SMatthew Knepley   }
6779566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&edgeMap));
6789566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
6797bee2925SMatthew Knepley   for (c = cStart; c < cEnd; ++c) {
6807bee2925SMatthew Knepley     PetscInt *closure = NULL;
6817bee2925SMatthew Knepley     PetscInt  clSize, cl, bval, fval;
6827bee2925SMatthew Knepley 
6839566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure));
6849566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(bodyLabel, c, &bval));
6859566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(faceLabel, c, &fval));
6867bee2925SMatthew Knepley     for (cl = 0; cl < clSize * 2; cl += 2) {
6879566063dSJacob Faibussowitsch       PetscCall(DMLabelSetValue(bodyLabel, closure[cl], bval));
6889566063dSJacob Faibussowitsch       PetscCall(DMLabelSetValue(faceLabel, closure[cl], fval));
6897bee2925SMatthew Knepley     }
6909566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure));
6917bee2925SMatthew Knepley   }
6927bee2925SMatthew Knepley   *newdm = dm;
6937bee2925SMatthew Knepley   PetscFunctionReturn(0);
6947bee2925SMatthew Knepley }
695c1cad2e7SMatthew G. Knepley 
696*9371c9d4SSatish Balay static PetscErrorCode DMPlexCreateEGADS(MPI_Comm comm, ego context, ego model, DM *newdm) {
697c1cad2e7SMatthew G. Knepley   DMLabel         bodyLabel, faceLabel, edgeLabel, vertexLabel;
698c1cad2e7SMatthew G. Knepley   // EGADS/EGADSLite variables
699c1cad2e7SMatthew G. Knepley   ego             geom, *bodies, *mobjs, *fobjs, *lobjs, *eobjs, *nobjs;
700c1cad2e7SMatthew G. Knepley   ego             topRef, prev, next;
701c1cad2e7SMatthew G. Knepley   int             oclass, mtype, nbodies, *senses, *lSenses, *eSenses;
702c1cad2e7SMatthew G. Knepley   int             b;
703c1cad2e7SMatthew G. Knepley   // PETSc variables
704c1cad2e7SMatthew G. Knepley   DM              dm;
705c1cad2e7SMatthew G. Knepley   PetscHMapI      edgeMap = NULL, bodyIndexMap = NULL, bodyVertexMap = NULL, bodyEdgeMap = NULL, bodyFaceMap = NULL, bodyEdgeGlobalMap = NULL;
706c1cad2e7SMatthew G. Knepley   PetscInt        dim = -1, cdim = -1, numCorners = 0, numVertices = 0, numEdges = 0, numFaces = 0, numCells = 0, edgeCntr = 0;
707c1cad2e7SMatthew G. Knepley   PetscInt        cellCntr = 0, numPoints = 0;
708c1cad2e7SMatthew G. Knepley   PetscInt       *cells  = NULL;
709c1cad2e7SMatthew G. Knepley   const PetscInt *cone   = NULL;
710c1cad2e7SMatthew G. Knepley   PetscReal      *coords = NULL;
711c1cad2e7SMatthew G. Knepley   PetscMPIInt     rank;
712c1cad2e7SMatthew G. Knepley 
713c1cad2e7SMatthew G. Knepley   PetscFunctionBeginUser;
7149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
715c5853193SPierre Jolivet   if (rank == 0) {
716c1cad2e7SMatthew G. Knepley     // ---------------------------------------------------------------------------------------------------
717c1cad2e7SMatthew G. Knepley     // Generate Petsc Plex
718c1cad2e7SMatthew G. Knepley     //  Get all Nodes in model, record coordinates in a correctly formatted array
719c1cad2e7SMatthew G. Knepley     //  Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array
720c1cad2e7SMatthew G. Knepley     //  We need to uniformly refine the initial geometry to guarantee a valid mesh
721c1cad2e7SMatthew G. Knepley 
722c1cad2e7SMatthew G. Knepley     // Caluculate cell and vertex sizes
7239566063dSJacob Faibussowitsch     PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses));
724c1cad2e7SMatthew G. Knepley 
7259566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&edgeMap));
7269566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&bodyIndexMap));
7279566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&bodyVertexMap));
7289566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&bodyEdgeMap));
7299566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&bodyEdgeGlobalMap));
7309566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&bodyFaceMap));
731c1cad2e7SMatthew G. Knepley 
732c1cad2e7SMatthew G. Knepley     for (b = 0; b < nbodies; ++b) {
733c1cad2e7SMatthew G. Knepley       ego           body = bodies[b];
734c1cad2e7SMatthew G. Knepley       int           Nf, Ne, Nv;
735c1cad2e7SMatthew G. Knepley       PetscHashIter BIiter, BViter, BEiter, BEGiter, BFiter, EMiter;
736c1cad2e7SMatthew G. Knepley       PetscBool     BIfound, BVfound, BEfound, BEGfound, BFfound, EMfound;
737c1cad2e7SMatthew G. Knepley 
7389566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(bodyIndexMap, b, &BIiter, &BIfound));
7399566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound));
7409566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound));
7419566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound));
7429566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound));
743c1cad2e7SMatthew G. Knepley 
7449566063dSJacob Faibussowitsch       if (!BIfound) PetscCall(PetscHMapISet(bodyIndexMap, b, numFaces + numEdges + numVertices));
7459566063dSJacob Faibussowitsch       if (!BVfound) PetscCall(PetscHMapISet(bodyVertexMap, b, numVertices));
7469566063dSJacob Faibussowitsch       if (!BEfound) PetscCall(PetscHMapISet(bodyEdgeMap, b, numEdges));
7479566063dSJacob Faibussowitsch       if (!BEGfound) PetscCall(PetscHMapISet(bodyEdgeGlobalMap, b, edgeCntr));
7489566063dSJacob Faibussowitsch       if (!BFfound) PetscCall(PetscHMapISet(bodyFaceMap, b, numFaces));
749c1cad2e7SMatthew G. Knepley 
7509566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
7519566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
7529566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
753c1cad2e7SMatthew G. Knepley       EG_free(fobjs);
754c1cad2e7SMatthew G. Knepley       EG_free(eobjs);
755c1cad2e7SMatthew G. Knepley       EG_free(nobjs);
756c1cad2e7SMatthew G. Knepley 
757c1cad2e7SMatthew G. Knepley       // Remove DEGENERATE EDGES from Edge count
7589566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
759c1cad2e7SMatthew G. Knepley       int Netemp = 0;
760c1cad2e7SMatthew G. Knepley       for (int e = 0; e < Ne; ++e) {
761c1cad2e7SMatthew G. Knepley         ego edge = eobjs[e];
762c1cad2e7SMatthew G. Knepley         int eid;
763c1cad2e7SMatthew G. Knepley 
7649566063dSJacob Faibussowitsch         PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
7655f80ce2aSJacob Faibussowitsch         eid = EG_indexBodyTopo(body, edge);
766c1cad2e7SMatthew G. Knepley 
7679566063dSJacob Faibussowitsch         PetscCall(PetscHMapIFind(edgeMap, edgeCntr + eid - 1, &EMiter, &EMfound));
768c1cad2e7SMatthew G. Knepley         if (mtype == DEGENERATE) {
7699566063dSJacob Faibussowitsch           if (!EMfound) PetscCall(PetscHMapISet(edgeMap, edgeCntr + eid - 1, -1));
770*9371c9d4SSatish Balay         } else {
771c1cad2e7SMatthew G. Knepley           ++Netemp;
7729566063dSJacob Faibussowitsch           if (!EMfound) PetscCall(PetscHMapISet(edgeMap, edgeCntr + eid - 1, Netemp));
773c1cad2e7SMatthew G. Knepley         }
774c1cad2e7SMatthew G. Knepley       }
775c1cad2e7SMatthew G. Knepley       EG_free(eobjs);
776c1cad2e7SMatthew G. Knepley 
777c1cad2e7SMatthew G. Knepley       // Determine Number of Cells
7789566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
779c1cad2e7SMatthew G. Knepley       for (int f = 0; f < Nf; ++f) {
780c1cad2e7SMatthew G. Knepley         ego face     = fobjs[f];
781c1cad2e7SMatthew G. Knepley         int edgeTemp = 0;
782c1cad2e7SMatthew G. Knepley 
7839566063dSJacob Faibussowitsch         PetscCall(EG_getBodyTopos(body, face, EDGE, &Ne, &eobjs));
784c1cad2e7SMatthew G. Knepley         for (int e = 0; e < Ne; ++e) {
785c1cad2e7SMatthew G. Knepley           ego edge = eobjs[e];
786c1cad2e7SMatthew G. Knepley 
7879566063dSJacob Faibussowitsch           PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
788c1cad2e7SMatthew G. Knepley           if (mtype != DEGENERATE) { ++edgeTemp; }
789c1cad2e7SMatthew G. Knepley         }
790c1cad2e7SMatthew G. Knepley         numCells += (2 * edgeTemp);
791c1cad2e7SMatthew G. Knepley         EG_free(eobjs);
792c1cad2e7SMatthew G. Knepley       }
793c1cad2e7SMatthew G. Knepley       EG_free(fobjs);
794c1cad2e7SMatthew G. Knepley 
795c1cad2e7SMatthew G. Knepley       numFaces += Nf;
796c1cad2e7SMatthew G. Knepley       numEdges += Netemp;
797c1cad2e7SMatthew G. Knepley       numVertices += Nv;
798c1cad2e7SMatthew G. Knepley       edgeCntr += Ne;
799c1cad2e7SMatthew G. Knepley     }
800c1cad2e7SMatthew G. Knepley 
801c1cad2e7SMatthew G. Knepley     // Set up basic DMPlex parameters
802c1cad2e7SMatthew G. Knepley     dim        = 2;                                 // Assumes 3D Models :: Need to handle 2D modles in the future
803c1cad2e7SMatthew G. Knepley     cdim       = 3;                                 // Assumes 3D Models :: Need to update to handle 2D modles in future
804c1cad2e7SMatthew G. Knepley     numCorners = 3;                                 // Split Faces into triangles
805c1cad2e7SMatthew G. Knepley     numPoints  = numVertices + numEdges + numFaces; // total number of coordinate points
806c1cad2e7SMatthew G. Knepley 
8079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(numPoints * cdim, &coords, numCells * numCorners, &cells));
808c1cad2e7SMatthew G. Knepley 
809c1cad2e7SMatthew G. Knepley     // Get Vertex Coordinates and Set up Cells
810c1cad2e7SMatthew G. Knepley     for (b = 0; b < nbodies; ++b) {
811c1cad2e7SMatthew G. Knepley       ego           body = bodies[b];
812c1cad2e7SMatthew G. Knepley       int           Nf, Ne, Nv;
813c1cad2e7SMatthew G. Knepley       PetscInt      bodyVertexIndexStart, bodyEdgeIndexStart, bodyEdgeGlobalIndexStart, bodyFaceIndexStart;
814c1cad2e7SMatthew G. Knepley       PetscHashIter BViter, BEiter, BEGiter, BFiter, EMiter;
815c1cad2e7SMatthew G. Knepley       PetscBool     BVfound, BEfound, BEGfound, BFfound, EMfound;
816c1cad2e7SMatthew G. Knepley 
817c1cad2e7SMatthew G. Knepley       // Vertices on Current Body
8189566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs));
819c1cad2e7SMatthew G. Knepley 
8209566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound));
82128b400f6SJacob Faibussowitsch       PetscCheck(BVfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyVertexMap", b);
8229566063dSJacob Faibussowitsch       PetscCall(PetscHMapIGet(bodyVertexMap, b, &bodyVertexIndexStart));
823c1cad2e7SMatthew G. Knepley 
824c1cad2e7SMatthew G. Knepley       for (int v = 0; v < Nv; ++v) {
825c1cad2e7SMatthew G. Knepley         ego    vertex = nobjs[v];
826c1cad2e7SMatthew G. Knepley         double limits[4];
827c1cad2e7SMatthew G. Knepley         int    id, dummy;
828c1cad2e7SMatthew G. Knepley 
8299566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses));
8305f80ce2aSJacob Faibussowitsch         id = EG_indexBodyTopo(body, vertex);
831c1cad2e7SMatthew G. Knepley 
832c1cad2e7SMatthew G. Knepley         coords[(bodyVertexIndexStart + id - 1) * cdim + 0] = limits[0];
833c1cad2e7SMatthew G. Knepley         coords[(bodyVertexIndexStart + id - 1) * cdim + 1] = limits[1];
834c1cad2e7SMatthew G. Knepley         coords[(bodyVertexIndexStart + id - 1) * cdim + 2] = limits[2];
835c1cad2e7SMatthew G. Knepley       }
836c1cad2e7SMatthew G. Knepley       EG_free(nobjs);
837c1cad2e7SMatthew G. Knepley 
838c1cad2e7SMatthew G. Knepley       // Edge Midpoint Vertices on Current Body
8399566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs));
840c1cad2e7SMatthew G. Knepley 
8419566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound));
84228b400f6SJacob Faibussowitsch       PetscCheck(BEfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeMap", b);
8439566063dSJacob Faibussowitsch       PetscCall(PetscHMapIGet(bodyEdgeMap, b, &bodyEdgeIndexStart));
844c1cad2e7SMatthew G. Knepley 
8459566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound));
84628b400f6SJacob Faibussowitsch       PetscCheck(BEGfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeGlobalMap", b);
8479566063dSJacob Faibussowitsch       PetscCall(PetscHMapIGet(bodyEdgeGlobalMap, b, &bodyEdgeGlobalIndexStart));
848c1cad2e7SMatthew G. Knepley 
849c1cad2e7SMatthew G. Knepley       for (int e = 0; e < Ne; ++e) {
850c1cad2e7SMatthew G. Knepley         ego    edge = eobjs[e];
851c1cad2e7SMatthew G. Knepley         double range[2], avgt[1], cntrPnt[9];
852c1cad2e7SMatthew G. Knepley         int    eid, eOffset;
853c1cad2e7SMatthew G. Knepley         int    periodic;
854c1cad2e7SMatthew G. Knepley 
8559566063dSJacob Faibussowitsch         PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
856c1cad2e7SMatthew G. Knepley         if (mtype == DEGENERATE) { continue; }
857c1cad2e7SMatthew G. Knepley 
8585f80ce2aSJacob Faibussowitsch         eid = EG_indexBodyTopo(body, edge);
859c1cad2e7SMatthew G. Knepley 
860c1cad2e7SMatthew G. Knepley         // get relative offset from globalEdgeID Vector
8619566063dSJacob Faibussowitsch         PetscCall(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound));
86228b400f6SJacob Faibussowitsch         PetscCheck(EMfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %d not found in edgeMap", bodyEdgeGlobalIndexStart + eid - 1);
8639566063dSJacob Faibussowitsch         PetscCall(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset));
864c1cad2e7SMatthew G. Knepley 
8659566063dSJacob Faibussowitsch         PetscCall(EG_getRange(edge, range, &periodic));
866c1cad2e7SMatthew G. Knepley         avgt[0] = (range[0] + range[1]) / 2.;
867c1cad2e7SMatthew G. Knepley 
8689566063dSJacob Faibussowitsch         PetscCall(EG_evaluate(edge, avgt, cntrPnt));
869c1cad2e7SMatthew G. Knepley         coords[(numVertices + bodyEdgeIndexStart + eOffset - 1) * cdim + 0] = cntrPnt[0];
870c1cad2e7SMatthew G. Knepley         coords[(numVertices + bodyEdgeIndexStart + eOffset - 1) * cdim + 1] = cntrPnt[1];
871c1cad2e7SMatthew G. Knepley         coords[(numVertices + bodyEdgeIndexStart + eOffset - 1) * cdim + 2] = cntrPnt[2];
872c1cad2e7SMatthew G. Knepley       }
873c1cad2e7SMatthew G. Knepley       EG_free(eobjs);
874c1cad2e7SMatthew G. Knepley 
875c1cad2e7SMatthew G. Knepley       // Face Midpoint Vertices on Current Body
8769566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
877c1cad2e7SMatthew G. Knepley 
8789566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound));
87928b400f6SJacob Faibussowitsch       PetscCheck(BFfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyFaceMap", b);
8809566063dSJacob Faibussowitsch       PetscCall(PetscHMapIGet(bodyFaceMap, b, &bodyFaceIndexStart));
881c1cad2e7SMatthew G. Knepley 
882c1cad2e7SMatthew G. Knepley       for (int f = 0; f < Nf; ++f) {
883c1cad2e7SMatthew G. Knepley         ego    face = fobjs[f];
884c1cad2e7SMatthew G. Knepley         double range[4], avgUV[2], cntrPnt[18];
885c1cad2e7SMatthew G. Knepley         int    peri, id;
886c1cad2e7SMatthew G. Knepley 
887c1cad2e7SMatthew G. Knepley         id = EG_indexBodyTopo(body, face);
8889566063dSJacob Faibussowitsch         PetscCall(EG_getRange(face, range, &peri));
889c1cad2e7SMatthew G. Knepley 
890c1cad2e7SMatthew G. Knepley         avgUV[0] = (range[0] + range[1]) / 2.;
891c1cad2e7SMatthew G. Knepley         avgUV[1] = (range[2] + range[3]) / 2.;
8929566063dSJacob Faibussowitsch         PetscCall(EG_evaluate(face, avgUV, cntrPnt));
893c1cad2e7SMatthew G. Knepley 
894c1cad2e7SMatthew G. Knepley         coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1) * cdim + 0] = cntrPnt[0];
895c1cad2e7SMatthew G. Knepley         coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1) * cdim + 1] = cntrPnt[1];
896c1cad2e7SMatthew G. Knepley         coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1) * cdim + 2] = cntrPnt[2];
897c1cad2e7SMatthew G. Knepley       }
898c1cad2e7SMatthew G. Knepley       EG_free(fobjs);
899c1cad2e7SMatthew G. Knepley 
900c1cad2e7SMatthew G. Knepley       // Define Cells :: Note - This could be incorporated in the Face Midpoint Vertices Loop but was kept separate for clarity
9019566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
902c1cad2e7SMatthew G. Knepley       for (int f = 0; f < Nf; ++f) {
903c1cad2e7SMatthew G. Knepley         ego face = fobjs[f];
904c1cad2e7SMatthew G. Knepley         int fID, midFaceID, midPntID, startID, endID, Nl;
905c1cad2e7SMatthew G. Knepley 
9065f80ce2aSJacob Faibussowitsch         fID       = EG_indexBodyTopo(body, face);
907c1cad2e7SMatthew G. Knepley         midFaceID = numVertices + numEdges + bodyFaceIndexStart + fID - 1;
908c1cad2e7SMatthew G. Knepley         // Must Traverse Loop to ensure we have all necessary information like the sense (+/- 1) of the edges.
909c1cad2e7SMatthew G. Knepley         // TODO :: Only handles single loop faces (No holes). The choices for handling multiloop faces are:
910c1cad2e7SMatthew G. Knepley         //            1) Use the DMPlexCreateEGADSFromFile() with the -dm_plex_egads_with_tess = 1 option.
911c1cad2e7SMatthew G. Knepley         //               This will use a default EGADS tessellation as an initial surface mesh.
912c1cad2e7SMatthew G. Knepley         //            2) Create the initial surface mesh via a 2D mesher :: Currently not availble (?future?)
913c1cad2e7SMatthew G. Knepley         //               May I suggest the XXXX as a starting point?
914c1cad2e7SMatthew G. Knepley 
9159566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(face, &geom, &oclass, &mtype, NULL, &Nl, &lobjs, &lSenses));
916c1cad2e7SMatthew G. Knepley 
91708401ef6SPierre Jolivet         PetscCheck(Nl <= 1, 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);
918c1cad2e7SMatthew G. Knepley         for (int l = 0; l < Nl; ++l) {
919c1cad2e7SMatthew G. Knepley           ego loop = lobjs[l];
920c1cad2e7SMatthew G. Knepley 
9219566063dSJacob Faibussowitsch           PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses));
922c1cad2e7SMatthew G. Knepley           for (int e = 0; e < Ne; ++e) {
923c1cad2e7SMatthew G. Knepley             ego edge = eobjs[e];
924c1cad2e7SMatthew G. Knepley             int eid, eOffset;
925c1cad2e7SMatthew G. Knepley 
9269566063dSJacob Faibussowitsch             PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
927c1cad2e7SMatthew G. Knepley             eid = EG_indexBodyTopo(body, edge);
928c1cad2e7SMatthew G. Knepley             if (mtype == DEGENERATE) { continue; }
929c1cad2e7SMatthew G. Knepley 
930c1cad2e7SMatthew G. Knepley             // get relative offset from globalEdgeID Vector
9319566063dSJacob Faibussowitsch             PetscCall(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound));
93228b400f6SJacob Faibussowitsch             PetscCheck(EMfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %d of Body %d not found in edgeMap. Global Edge ID :: %d", eid, b, bodyEdgeGlobalIndexStart + eid - 1);
9339566063dSJacob Faibussowitsch             PetscCall(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset));
934c1cad2e7SMatthew G. Knepley 
935c1cad2e7SMatthew G. Knepley             midPntID = numVertices + bodyEdgeIndexStart + eOffset - 1;
936c1cad2e7SMatthew G. Knepley 
9379566063dSJacob Faibussowitsch             PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses));
938c1cad2e7SMatthew G. Knepley 
939*9371c9d4SSatish Balay             if (eSenses[e] > 0) {
940*9371c9d4SSatish Balay               startID = EG_indexBodyTopo(body, nobjs[0]);
941*9371c9d4SSatish Balay               endID   = EG_indexBodyTopo(body, nobjs[1]);
942*9371c9d4SSatish Balay             } else {
943*9371c9d4SSatish Balay               startID = EG_indexBodyTopo(body, nobjs[1]);
944*9371c9d4SSatish Balay               endID   = EG_indexBodyTopo(body, nobjs[0]);
945*9371c9d4SSatish Balay             }
946c1cad2e7SMatthew G. Knepley 
947c1cad2e7SMatthew G. Knepley             // Define 2 Cells per Edge with correct orientation
948c1cad2e7SMatthew G. Knepley             cells[cellCntr * numCorners + 0] = midFaceID;
949c1cad2e7SMatthew G. Knepley             cells[cellCntr * numCorners + 1] = bodyVertexIndexStart + startID - 1;
950c1cad2e7SMatthew G. Knepley             cells[cellCntr * numCorners + 2] = midPntID;
951c1cad2e7SMatthew G. Knepley 
952c1cad2e7SMatthew G. Knepley             cells[cellCntr * numCorners + 3] = midFaceID;
953c1cad2e7SMatthew G. Knepley             cells[cellCntr * numCorners + 4] = midPntID;
954c1cad2e7SMatthew G. Knepley             cells[cellCntr * numCorners + 5] = bodyVertexIndexStart + endID - 1;
955c1cad2e7SMatthew G. Knepley 
956c1cad2e7SMatthew G. Knepley             cellCntr = cellCntr + 2;
957c1cad2e7SMatthew G. Knepley           }
958c1cad2e7SMatthew G. Knepley         }
959c1cad2e7SMatthew G. Knepley       }
960c1cad2e7SMatthew G. Knepley       EG_free(fobjs);
961c1cad2e7SMatthew G. Knepley     }
962c1cad2e7SMatthew G. Knepley   }
963c1cad2e7SMatthew G. Knepley 
964c1cad2e7SMatthew G. Knepley   // Generate DMPlex
9659566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numPoints, numCorners, PETSC_TRUE, cells, cdim, coords, &dm));
9669566063dSJacob Faibussowitsch   PetscCall(PetscFree2(coords, cells));
96763a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(dm, " Total Number of Unique Cells    = %" PetscInt_FMT " \n", numCells));
96863a3b9bcSJacob Faibussowitsch   PetscCall(PetscInfo(dm, " Total Number of Unique Vertices = %" PetscInt_FMT " \n", numVertices));
969c1cad2e7SMatthew G. Knepley 
970c1cad2e7SMatthew G. Knepley   // Embed EGADS model in DM
971c1cad2e7SMatthew G. Knepley   {
972c1cad2e7SMatthew G. Knepley     PetscContainer modelObj, contextObj;
973c1cad2e7SMatthew G. Knepley 
9749566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &modelObj));
9759566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(modelObj, model));
9769566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Model", (PetscObject)modelObj));
9779566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&modelObj));
978c1cad2e7SMatthew G. Knepley 
9799566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &contextObj));
9809566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(contextObj, context));
9819566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private));
9829566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Context", (PetscObject)contextObj));
9839566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&contextObj));
984c1cad2e7SMatthew G. Knepley   }
985c1cad2e7SMatthew G. Knepley   // Label points
986c1cad2e7SMatthew G. Knepley   PetscInt nStart, nEnd;
987c1cad2e7SMatthew G. Knepley 
9889566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Body ID"));
9899566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
9909566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Face ID"));
9919566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
9929566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Edge ID"));
9939566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
9949566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Vertex ID"));
9959566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
996c1cad2e7SMatthew G. Knepley 
9979566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 2, &nStart, &nEnd));
998c1cad2e7SMatthew G. Knepley 
999c1cad2e7SMatthew G. Knepley   cellCntr = 0;
1000c1cad2e7SMatthew G. Knepley   for (b = 0; b < nbodies; ++b) {
1001c1cad2e7SMatthew G. Knepley     ego           body = bodies[b];
1002c1cad2e7SMatthew G. Knepley     int           Nv, Ne, Nf;
1003c1cad2e7SMatthew G. Knepley     PetscInt      bodyVertexIndexStart, bodyEdgeIndexStart, bodyEdgeGlobalIndexStart, bodyFaceIndexStart;
1004c1cad2e7SMatthew G. Knepley     PetscHashIter BViter, BEiter, BEGiter, BFiter, EMiter;
1005c1cad2e7SMatthew G. Knepley     PetscBool     BVfound, BEfound, BEGfound, BFfound, EMfound;
1006c1cad2e7SMatthew G. Knepley 
10079566063dSJacob Faibussowitsch     PetscCall(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound));
100828b400f6SJacob Faibussowitsch     PetscCheck(BVfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyVertexMap", b);
10099566063dSJacob Faibussowitsch     PetscCall(PetscHMapIGet(bodyVertexMap, b, &bodyVertexIndexStart));
1010c1cad2e7SMatthew G. Knepley 
10119566063dSJacob Faibussowitsch     PetscCall(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound));
101228b400f6SJacob Faibussowitsch     PetscCheck(BEfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeMap", b);
10139566063dSJacob Faibussowitsch     PetscCall(PetscHMapIGet(bodyEdgeMap, b, &bodyEdgeIndexStart));
1014c1cad2e7SMatthew G. Knepley 
10159566063dSJacob Faibussowitsch     PetscCall(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound));
101628b400f6SJacob Faibussowitsch     PetscCheck(BFfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyFaceMap", b);
10179566063dSJacob Faibussowitsch     PetscCall(PetscHMapIGet(bodyFaceMap, b, &bodyFaceIndexStart));
1018c1cad2e7SMatthew G. Knepley 
10199566063dSJacob Faibussowitsch     PetscCall(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound));
102028b400f6SJacob Faibussowitsch     PetscCheck(BEGfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeGlobalMap", b);
10219566063dSJacob Faibussowitsch     PetscCall(PetscHMapIGet(bodyEdgeGlobalMap, b, &bodyEdgeGlobalIndexStart));
1022c1cad2e7SMatthew G. Knepley 
10239566063dSJacob Faibussowitsch     PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1024c1cad2e7SMatthew G. Knepley     for (int f = 0; f < Nf; ++f) {
1025c1cad2e7SMatthew G. Knepley       ego face = fobjs[f];
1026c1cad2e7SMatthew G. Knepley       int fID, Nl;
1027c1cad2e7SMatthew G. Knepley 
10285f80ce2aSJacob Faibussowitsch       fID = EG_indexBodyTopo(body, face);
1029c1cad2e7SMatthew G. Knepley 
10309566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, face, LOOP, &Nl, &lobjs));
1031c1cad2e7SMatthew G. Knepley       for (int l = 0; l < Nl; ++l) {
1032c1cad2e7SMatthew G. Knepley         ego loop = lobjs[l];
1033c1cad2e7SMatthew G. Knepley         int lid;
1034c1cad2e7SMatthew G. Knepley 
10355f80ce2aSJacob Faibussowitsch         lid = EG_indexBodyTopo(body, loop);
103608401ef6SPierre Jolivet         PetscCheck(Nl <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf);
1037c1cad2e7SMatthew G. Knepley 
10389566063dSJacob Faibussowitsch         PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses));
1039c1cad2e7SMatthew G. Knepley         for (int e = 0; e < Ne; ++e) {
1040c1cad2e7SMatthew G. Knepley           ego edge = eobjs[e];
1041c1cad2e7SMatthew G. Knepley           int eid, eOffset;
1042c1cad2e7SMatthew G. Knepley 
1043c1cad2e7SMatthew G. Knepley           // Skip DEGENERATE Edges
10449566063dSJacob Faibussowitsch           PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next));
1045c1cad2e7SMatthew G. Knepley           if (mtype == DEGENERATE) { continue; }
10465f80ce2aSJacob Faibussowitsch           eid = EG_indexBodyTopo(body, edge);
1047c1cad2e7SMatthew G. Knepley 
1048c1cad2e7SMatthew G. Knepley           // get relative offset from globalEdgeID Vector
10499566063dSJacob Faibussowitsch           PetscCall(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound));
105028b400f6SJacob Faibussowitsch           PetscCheck(EMfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %d of Body %d not found in edgeMap. Global Edge ID :: %d", eid, b, bodyEdgeGlobalIndexStart + eid - 1);
10519566063dSJacob Faibussowitsch           PetscCall(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset));
1052c1cad2e7SMatthew G. Knepley 
10539566063dSJacob Faibussowitsch           PetscCall(EG_getBodyTopos(body, edge, NODE, &Nv, &nobjs));
1054c1cad2e7SMatthew G. Knepley           for (int v = 0; v < Nv; ++v) {
1055c1cad2e7SMatthew G. Knepley             ego vertex = nobjs[v];
1056c1cad2e7SMatthew G. Knepley             int vID;
1057c1cad2e7SMatthew G. Knepley 
10585f80ce2aSJacob Faibussowitsch             vID = EG_indexBodyTopo(body, vertex);
10599566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(bodyLabel, nStart + bodyVertexIndexStart + vID - 1, b));
10609566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(vertexLabel, nStart + bodyVertexIndexStart + vID - 1, vID));
1061c1cad2e7SMatthew G. Knepley           }
1062c1cad2e7SMatthew G. Knepley           EG_free(nobjs);
1063c1cad2e7SMatthew G. Knepley 
10649566063dSJacob Faibussowitsch           PetscCall(DMLabelSetValue(bodyLabel, nStart + numVertices + bodyEdgeIndexStart + eOffset - 1, b));
10659566063dSJacob Faibussowitsch           PetscCall(DMLabelSetValue(edgeLabel, nStart + numVertices + bodyEdgeIndexStart + eOffset - 1, eid));
1066c1cad2e7SMatthew G. Knepley 
1067c1cad2e7SMatthew G. Knepley           // Define Cell faces
1068c1cad2e7SMatthew G. Knepley           for (int jj = 0; jj < 2; ++jj) {
10699566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(bodyLabel, cellCntr, b));
10709566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(faceLabel, cellCntr, fID));
10719566063dSJacob Faibussowitsch             PetscCall(DMPlexGetCone(dm, cellCntr, &cone));
1072c1cad2e7SMatthew G. Knepley 
10739566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(bodyLabel, cone[0], b));
10749566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(faceLabel, cone[0], fID));
1075c1cad2e7SMatthew G. Knepley 
10769566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(bodyLabel, cone[1], b));
10779566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(edgeLabel, cone[1], eid));
1078c1cad2e7SMatthew G. Knepley 
10799566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(bodyLabel, cone[2], b));
10809566063dSJacob Faibussowitsch             PetscCall(DMLabelSetValue(faceLabel, cone[2], fID));
1081c1cad2e7SMatthew G. Knepley 
1082c1cad2e7SMatthew G. Knepley             cellCntr = cellCntr + 1;
1083c1cad2e7SMatthew G. Knepley           }
1084c1cad2e7SMatthew G. Knepley         }
1085c1cad2e7SMatthew G. Knepley       }
1086c1cad2e7SMatthew G. Knepley       EG_free(lobjs);
1087c1cad2e7SMatthew G. Knepley 
10889566063dSJacob Faibussowitsch       PetscCall(DMLabelSetValue(bodyLabel, nStart + numVertices + numEdges + bodyFaceIndexStart + fID - 1, b));
10899566063dSJacob Faibussowitsch       PetscCall(DMLabelSetValue(faceLabel, nStart + numVertices + numEdges + bodyFaceIndexStart + fID - 1, fID));
1090c1cad2e7SMatthew G. Knepley     }
1091c1cad2e7SMatthew G. Knepley     EG_free(fobjs);
1092c1cad2e7SMatthew G. Knepley   }
1093c1cad2e7SMatthew G. Knepley 
10949566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&edgeMap));
10959566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&bodyIndexMap));
10969566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&bodyVertexMap));
10979566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&bodyEdgeMap));
10989566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&bodyEdgeGlobalMap));
10999566063dSJacob Faibussowitsch   PetscCall(PetscHMapIDestroy(&bodyFaceMap));
1100c1cad2e7SMatthew G. Knepley 
1101c1cad2e7SMatthew G. Knepley   *newdm = dm;
1102c1cad2e7SMatthew G. Knepley   PetscFunctionReturn(0);
1103c1cad2e7SMatthew G. Knepley }
1104c1cad2e7SMatthew G. Knepley 
1105*9371c9d4SSatish Balay static PetscErrorCode DMPlexCreateEGADS_Tess_Internal(MPI_Comm comm, ego context, ego model, DM *newdm) {
1106c1cad2e7SMatthew G. Knepley   DMLabel         bodyLabel, faceLabel, edgeLabel, vertexLabel;
1107c1cad2e7SMatthew G. Knepley   /* EGADSLite variables */
1108c1cad2e7SMatthew G. Knepley   ego             geom, *bodies, *fobjs;
1109c1cad2e7SMatthew G. Knepley   int             b, oclass, mtype, nbodies, *senses;
1110c1cad2e7SMatthew G. Knepley   int             totalNumTris = 0, totalNumPoints = 0;
1111c1cad2e7SMatthew G. Knepley   double          boundBox[6] = {0., 0., 0., 0., 0., 0.}, tessSize;
1112c1cad2e7SMatthew G. Knepley   /* PETSc variables */
1113c1cad2e7SMatthew G. Knepley   DM              dm;
1114c1cad2e7SMatthew G. Knepley   PetscHMapI      pointIndexStartMap = NULL, triIndexStartMap = NULL, pTypeLabelMap = NULL, pIndexLabelMap = NULL;
1115c1cad2e7SMatthew G. Knepley   PetscHMapI      pBodyIndexLabelMap = NULL, triFaceIDLabelMap = NULL, triBodyIDLabelMap = NULL;
1116c1cad2e7SMatthew G. Knepley   PetscInt        dim = -1, cdim = -1, numCorners = 0, counter = 0;
1117c1cad2e7SMatthew G. Knepley   PetscInt       *cells  = NULL;
1118c1cad2e7SMatthew G. Knepley   const PetscInt *cone   = NULL;
1119c1cad2e7SMatthew G. Knepley   PetscReal      *coords = NULL;
1120c1cad2e7SMatthew G. Knepley   PetscMPIInt     rank;
1121c1cad2e7SMatthew G. Knepley 
1122c1cad2e7SMatthew G. Knepley   PetscFunctionBeginUser;
11239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1124c5853193SPierre Jolivet   if (rank == 0) {
1125c1cad2e7SMatthew G. Knepley     // ---------------------------------------------------------------------------------------------------
1126c1cad2e7SMatthew G. Knepley     // Generate Petsc Plex from EGADSlite created Tessellation of geometry
1127c1cad2e7SMatthew G. Knepley     // ---------------------------------------------------------------------------------------------------
1128c1cad2e7SMatthew G. Knepley 
1129c1cad2e7SMatthew G. Knepley     // Caluculate cell and vertex sizes
11309566063dSJacob Faibussowitsch     PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses));
1131c1cad2e7SMatthew G. Knepley 
11329566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&pointIndexStartMap));
11339566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&triIndexStartMap));
11349566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&pTypeLabelMap));
11359566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&pIndexLabelMap));
11369566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&pBodyIndexLabelMap));
11379566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&triFaceIDLabelMap));
11389566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&triBodyIDLabelMap));
1139c1cad2e7SMatthew G. Knepley 
1140c1cad2e7SMatthew G. Knepley     /* Create Tessellation of Bodies */
1141c1cad2e7SMatthew G. Knepley     ego tessArray[nbodies];
1142c1cad2e7SMatthew G. Knepley 
1143c1cad2e7SMatthew G. Knepley     for (b = 0; b < nbodies; ++b) {
1144c1cad2e7SMatthew G. Knepley       ego           body      = bodies[b];
1145c1cad2e7SMatthew G. Knepley       double        params[3] = {0.0, 0.0, 0.0}; // Parameters for Tessellation
1146c1cad2e7SMatthew G. Knepley       int           Nf, bodyNumPoints = 0, bodyNumTris = 0;
1147c1cad2e7SMatthew G. Knepley       PetscHashIter PISiter, TISiter;
1148c1cad2e7SMatthew G. Knepley       PetscBool     PISfound, TISfound;
1149c1cad2e7SMatthew G. Knepley 
1150c1cad2e7SMatthew G. Knepley       /* Store Start Index for each Body's Point and Tris */
11519566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(pointIndexStartMap, b, &PISiter, &PISfound));
11529566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(triIndexStartMap, b, &TISiter, &TISfound));
1153c1cad2e7SMatthew G. Knepley 
11549566063dSJacob Faibussowitsch       if (!PISfound) PetscCall(PetscHMapISet(pointIndexStartMap, b, totalNumPoints));
11559566063dSJacob Faibussowitsch       if (!TISfound) PetscCall(PetscHMapISet(triIndexStartMap, b, totalNumTris));
1156c1cad2e7SMatthew G. Knepley 
1157c1cad2e7SMatthew G. Knepley       /* Calculate Tessellation parameters based on Bounding Box */
1158c1cad2e7SMatthew G. Knepley       /* Get Bounding Box Dimensions of the BODY */
11599566063dSJacob Faibussowitsch       PetscCall(EG_getBoundingBox(body, boundBox));
1160c1cad2e7SMatthew G. Knepley       tessSize = boundBox[3] - boundBox[0];
1161c1cad2e7SMatthew G. Knepley       if (tessSize < boundBox[4] - boundBox[1]) tessSize = boundBox[4] - boundBox[1];
1162c1cad2e7SMatthew G. Knepley       if (tessSize < boundBox[5] - boundBox[2]) tessSize = boundBox[5] - boundBox[2];
1163c1cad2e7SMatthew G. Knepley 
1164c1cad2e7SMatthew G. Knepley       // TODO :: May want to give users tessellation parameter options //
1165c1cad2e7SMatthew G. Knepley       params[0] = 0.0250 * tessSize;
1166c1cad2e7SMatthew G. Knepley       params[1] = 0.0075 * tessSize;
1167c1cad2e7SMatthew G. Knepley       params[2] = 15.0;
1168c1cad2e7SMatthew G. Knepley 
11699566063dSJacob Faibussowitsch       PetscCall(EG_makeTessBody(body, params, &tessArray[b]));
1170c1cad2e7SMatthew G. Knepley 
11719566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1172c1cad2e7SMatthew G. Knepley 
1173c1cad2e7SMatthew G. Knepley       for (int f = 0; f < Nf; ++f) {
1174c1cad2e7SMatthew G. Knepley         ego           face = fobjs[f];
1175c1cad2e7SMatthew G. Knepley         int           len, fID, ntris;
1176c1cad2e7SMatthew G. Knepley         const int    *ptype, *pindex, *ptris, *ptric;
1177c1cad2e7SMatthew G. Knepley         const double *pxyz, *puv;
1178c1cad2e7SMatthew G. Knepley 
1179c1cad2e7SMatthew G. Knepley         // Get Face ID //
1180c1cad2e7SMatthew G. Knepley         fID = EG_indexBodyTopo(body, face);
1181c1cad2e7SMatthew G. Knepley 
1182c1cad2e7SMatthew G. Knepley         // Checkout the Surface Tessellation //
11839566063dSJacob Faibussowitsch         PetscCall(EG_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric));
1184c1cad2e7SMatthew G. Knepley 
1185c1cad2e7SMatthew G. Knepley         // Determine total number of triangle cells in the tessellation //
1186c1cad2e7SMatthew G. Knepley         bodyNumTris += (int)ntris;
1187c1cad2e7SMatthew G. Knepley 
1188c1cad2e7SMatthew G. Knepley         // Check out the point index and coordinate //
1189c1cad2e7SMatthew G. Knepley         for (int p = 0; p < len; ++p) {
1190c1cad2e7SMatthew G. Knepley           int global;
1191c1cad2e7SMatthew G. Knepley 
11929566063dSJacob Faibussowitsch           PetscCall(EG_localToGlobal(tessArray[b], fID, p + 1, &global));
1193c1cad2e7SMatthew G. Knepley 
1194c1cad2e7SMatthew G. Knepley           // Determine the total number of points in the tessellation //
1195c1cad2e7SMatthew G. Knepley           bodyNumPoints = PetscMax(bodyNumPoints, global);
1196c1cad2e7SMatthew G. Knepley         }
1197c1cad2e7SMatthew G. Knepley       }
1198c1cad2e7SMatthew G. Knepley       EG_free(fobjs);
1199c1cad2e7SMatthew G. Knepley 
1200c1cad2e7SMatthew G. Knepley       totalNumPoints += bodyNumPoints;
1201c1cad2e7SMatthew G. Knepley       totalNumTris += bodyNumTris;
1202c1cad2e7SMatthew G. Knepley     }
1203c5853193SPierre Jolivet     //}  - Original End of (rank == 0)
1204c1cad2e7SMatthew G. Knepley 
1205c1cad2e7SMatthew G. Knepley     dim        = 2;
1206c1cad2e7SMatthew G. Knepley     cdim       = 3;
1207c1cad2e7SMatthew G. Knepley     numCorners = 3;
1208c1cad2e7SMatthew G. Knepley     //PetscInt counter = 0;
1209c1cad2e7SMatthew G. Knepley 
1210c1cad2e7SMatthew G. Knepley     /* NEED TO DEFINE MATRICES/VECTORS TO STORE GEOM REFERENCE DATA   */
1211c1cad2e7SMatthew G. Knepley     /* Fill in below and use to define DMLabels after DMPlex creation */
12129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(totalNumPoints * cdim, &coords, totalNumTris * numCorners, &cells));
1213c1cad2e7SMatthew G. Knepley 
1214c1cad2e7SMatthew G. Knepley     for (b = 0; b < nbodies; ++b) {
1215c1cad2e7SMatthew G. Knepley       ego           body = bodies[b];
1216c1cad2e7SMatthew G. Knepley       int           Nf;
1217c1cad2e7SMatthew G. Knepley       PetscInt      pointIndexStart;
1218c1cad2e7SMatthew G. Knepley       PetscHashIter PISiter;
1219c1cad2e7SMatthew G. Knepley       PetscBool     PISfound;
1220c1cad2e7SMatthew G. Knepley 
12219566063dSJacob Faibussowitsch       PetscCall(PetscHMapIFind(pointIndexStartMap, b, &PISiter, &PISfound));
122228b400f6SJacob Faibussowitsch       PetscCheck(PISfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in pointIndexStartMap", b);
12239566063dSJacob Faibussowitsch       PetscCall(PetscHMapIGet(pointIndexStartMap, b, &pointIndexStart));
1224c1cad2e7SMatthew G. Knepley 
12259566063dSJacob Faibussowitsch       PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs));
1226c1cad2e7SMatthew G. Knepley 
1227c1cad2e7SMatthew G. Knepley       for (int f = 0; f < Nf; ++f) {
1228c1cad2e7SMatthew G. Knepley         /* Get Face Object */
1229c1cad2e7SMatthew G. Knepley         ego           face = fobjs[f];
1230c1cad2e7SMatthew G. Knepley         int           len, fID, ntris;
1231c1cad2e7SMatthew G. Knepley         const int    *ptype, *pindex, *ptris, *ptric;
1232c1cad2e7SMatthew G. Knepley         const double *pxyz, *puv;
1233c1cad2e7SMatthew G. Knepley 
1234c1cad2e7SMatthew G. Knepley         /* Get Face ID */
1235c1cad2e7SMatthew G. Knepley         fID = EG_indexBodyTopo(body, face);
1236c1cad2e7SMatthew G. Knepley 
1237c1cad2e7SMatthew G. Knepley         /* Checkout the Surface Tessellation */
12389566063dSJacob Faibussowitsch         PetscCall(EG_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric));
1239c1cad2e7SMatthew G. Knepley 
1240c1cad2e7SMatthew G. Knepley         /* Check out the point index and coordinate */
1241c1cad2e7SMatthew G. Knepley         for (int p = 0; p < len; ++p) {
1242c1cad2e7SMatthew G. Knepley           int           global;
1243c1cad2e7SMatthew G. Knepley           PetscHashIter PTLiter, PILiter, PBLiter;
1244c1cad2e7SMatthew G. Knepley           PetscBool     PTLfound, PILfound, PBLfound;
1245c1cad2e7SMatthew G. Knepley 
12469566063dSJacob Faibussowitsch           PetscCall(EG_localToGlobal(tessArray[b], fID, p + 1, &global));
1247c1cad2e7SMatthew G. Knepley 
1248c1cad2e7SMatthew G. Knepley           /* Set the coordinates array for DAG */
1249c1cad2e7SMatthew G. Knepley           coords[((global - 1 + pointIndexStart) * 3) + 0] = pxyz[(p * 3) + 0];
1250c1cad2e7SMatthew G. Knepley           coords[((global - 1 + pointIndexStart) * 3) + 1] = pxyz[(p * 3) + 1];
1251c1cad2e7SMatthew G. Knepley           coords[((global - 1 + pointIndexStart) * 3) + 2] = pxyz[(p * 3) + 2];
1252c1cad2e7SMatthew G. Knepley 
1253c1cad2e7SMatthew G. Knepley           /* Store Geometry Label Information for DMLabel assignment later */
12549566063dSJacob Faibussowitsch           PetscCall(PetscHMapIFind(pTypeLabelMap, global - 1 + pointIndexStart, &PTLiter, &PTLfound));
12559566063dSJacob Faibussowitsch           PetscCall(PetscHMapIFind(pIndexLabelMap, global - 1 + pointIndexStart, &PILiter, &PILfound));
12569566063dSJacob Faibussowitsch           PetscCall(PetscHMapIFind(pBodyIndexLabelMap, global - 1 + pointIndexStart, &PBLiter, &PBLfound));
1257c1cad2e7SMatthew G. Knepley 
12589566063dSJacob Faibussowitsch           if (!PTLfound) PetscCall(PetscHMapISet(pTypeLabelMap, global - 1 + pointIndexStart, ptype[p]));
12599566063dSJacob Faibussowitsch           if (!PILfound) PetscCall(PetscHMapISet(pIndexLabelMap, global - 1 + pointIndexStart, pindex[p]));
12609566063dSJacob Faibussowitsch           if (!PBLfound) PetscCall(PetscHMapISet(pBodyIndexLabelMap, global - 1 + pointIndexStart, b));
1261c1cad2e7SMatthew G. Knepley 
12629566063dSJacob Faibussowitsch           if (ptype[p] < 0) PetscCall(PetscHMapISet(pIndexLabelMap, global - 1 + pointIndexStart, fID));
1263c1cad2e7SMatthew G. Knepley         }
1264c1cad2e7SMatthew G. Knepley 
1265c1cad2e7SMatthew G. Knepley         for (int t = 0; t < (int)ntris; ++t) {
1266c1cad2e7SMatthew G. Knepley           int           global, globalA, globalB;
1267c1cad2e7SMatthew G. Knepley           PetscHashIter TFLiter, TBLiter;
1268c1cad2e7SMatthew G. Knepley           PetscBool     TFLfound, TBLfound;
1269c1cad2e7SMatthew G. Knepley 
12709566063dSJacob Faibussowitsch           PetscCall(EG_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 0], &global));
1271c1cad2e7SMatthew G. Knepley           cells[(counter * 3) + 0] = global - 1 + pointIndexStart;
1272c1cad2e7SMatthew G. Knepley 
12739566063dSJacob Faibussowitsch           PetscCall(EG_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 1], &globalA));
1274c1cad2e7SMatthew G. Knepley           cells[(counter * 3) + 1] = globalA - 1 + pointIndexStart;
1275c1cad2e7SMatthew G. Knepley 
12769566063dSJacob Faibussowitsch           PetscCall(EG_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 2], &globalB));
1277c1cad2e7SMatthew G. Knepley           cells[(counter * 3) + 2] = globalB - 1 + pointIndexStart;
1278c1cad2e7SMatthew G. Knepley 
12799566063dSJacob Faibussowitsch           PetscCall(PetscHMapIFind(triFaceIDLabelMap, counter, &TFLiter, &TFLfound));
12809566063dSJacob Faibussowitsch           PetscCall(PetscHMapIFind(triBodyIDLabelMap, counter, &TBLiter, &TBLfound));
1281c1cad2e7SMatthew G. Knepley 
12829566063dSJacob Faibussowitsch           if (!TFLfound) PetscCall(PetscHMapISet(triFaceIDLabelMap, counter, fID));
12839566063dSJacob Faibussowitsch           if (!TBLfound) PetscCall(PetscHMapISet(triBodyIDLabelMap, counter, b));
1284c1cad2e7SMatthew G. Knepley 
1285c1cad2e7SMatthew G. Knepley           counter += 1;
1286c1cad2e7SMatthew G. Knepley         }
1287c1cad2e7SMatthew G. Knepley       }
1288c1cad2e7SMatthew G. Knepley       EG_free(fobjs);
1289c1cad2e7SMatthew G. Knepley     }
1290c1cad2e7SMatthew G. Knepley   }
1291c1cad2e7SMatthew G. Knepley 
1292c1cad2e7SMatthew G. Knepley   //Build DMPlex
12939566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, totalNumTris, totalNumPoints, numCorners, PETSC_TRUE, cells, cdim, coords, &dm));
12949566063dSJacob Faibussowitsch   PetscCall(PetscFree2(coords, cells));
1295c1cad2e7SMatthew G. Knepley 
1296c1cad2e7SMatthew G. Knepley   // Embed EGADS model in DM
1297c1cad2e7SMatthew G. Knepley   {
1298c1cad2e7SMatthew G. Knepley     PetscContainer modelObj, contextObj;
1299c1cad2e7SMatthew G. Knepley 
13009566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &modelObj));
13019566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(modelObj, model));
13029566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Model", (PetscObject)modelObj));
13039566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&modelObj));
1304c1cad2e7SMatthew G. Knepley 
13059566063dSJacob Faibussowitsch     PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &contextObj));
13069566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetPointer(contextObj, context));
13079566063dSJacob Faibussowitsch     PetscCall(PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private));
13089566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Context", (PetscObject)contextObj));
13099566063dSJacob Faibussowitsch     PetscCall(PetscContainerDestroy(&contextObj));
1310c1cad2e7SMatthew G. Knepley   }
1311c1cad2e7SMatthew G. Knepley 
1312c1cad2e7SMatthew G. Knepley   // Label Points
13139566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Body ID"));
13149566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
13159566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Face ID"));
13169566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
13179566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Edge ID"));
13189566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
13199566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(dm, "EGADS Vertex ID"));
13209566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
1321c1cad2e7SMatthew G. Knepley 
1322c1cad2e7SMatthew G. Knepley   /* Get Number of DAG Nodes at each level */
1323c1cad2e7SMatthew G. Knepley   int fStart, fEnd, eStart, eEnd, nStart, nEnd;
1324c1cad2e7SMatthew G. Knepley 
13259566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &fStart, &fEnd));
13269566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 1, &eStart, &eEnd));
13279566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 2, &nStart, &nEnd));
1328c1cad2e7SMatthew G. Knepley 
1329c1cad2e7SMatthew G. Knepley   /* Set DMLabels for NODES */
1330c1cad2e7SMatthew G. Knepley   for (int n = nStart; n < nEnd; ++n) {
1331c1cad2e7SMatthew G. Knepley     int           pTypeVal, pIndexVal, pBodyVal;
1332c1cad2e7SMatthew G. Knepley     PetscHashIter PTLiter, PILiter, PBLiter;
1333c1cad2e7SMatthew G. Knepley     PetscBool     PTLfound, PILfound, PBLfound;
1334c1cad2e7SMatthew G. Knepley 
1335c1cad2e7SMatthew G. Knepley     //Converted to Hash Tables
13369566063dSJacob Faibussowitsch     PetscCall(PetscHMapIFind(pTypeLabelMap, n - nStart, &PTLiter, &PTLfound));
133728b400f6SJacob Faibussowitsch     PetscCheck(PTLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pTypeLabelMap", n);
13389566063dSJacob Faibussowitsch     PetscCall(PetscHMapIGet(pTypeLabelMap, n - nStart, &pTypeVal));
1339c1cad2e7SMatthew G. Knepley 
13409566063dSJacob Faibussowitsch     PetscCall(PetscHMapIFind(pIndexLabelMap, n - nStart, &PILiter, &PILfound));
134128b400f6SJacob Faibussowitsch     PetscCheck(PILfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pIndexLabelMap", n);
13429566063dSJacob Faibussowitsch     PetscCall(PetscHMapIGet(pIndexLabelMap, n - nStart, &pIndexVal));
1343c1cad2e7SMatthew G. Knepley 
13449566063dSJacob Faibussowitsch     PetscCall(PetscHMapIFind(pBodyIndexLabelMap, n - nStart, &PBLiter, &PBLfound));
134528b400f6SJacob Faibussowitsch     PetscCheck(PBLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pBodyLabelMap", n);
13469566063dSJacob Faibussowitsch     PetscCall(PetscHMapIGet(pBodyIndexLabelMap, n - nStart, &pBodyVal));
1347c1cad2e7SMatthew G. Knepley 
13489566063dSJacob Faibussowitsch     PetscCall(DMLabelSetValue(bodyLabel, n, pBodyVal));
13499566063dSJacob Faibussowitsch     if (pTypeVal == 0) PetscCall(DMLabelSetValue(vertexLabel, n, pIndexVal));
13509566063dSJacob Faibussowitsch     if (pTypeVal > 0) PetscCall(DMLabelSetValue(edgeLabel, n, pIndexVal));
13519566063dSJacob Faibussowitsch     if (pTypeVal < 0) PetscCall(DMLabelSetValue(faceLabel, n, pIndexVal));
1352c1cad2e7SMatthew G. Knepley   }
1353c1cad2e7SMatthew G. Knepley 
1354c1cad2e7SMatthew G. Knepley   /* Set DMLabels for Edges - Based on the DMLabels of the EDGE's NODES */
1355c1cad2e7SMatthew G. Knepley   for (int e = eStart; e < eEnd; ++e) {
1356c1cad2e7SMatthew G. Knepley     int bodyID_0, vertexID_0, vertexID_1, edgeID_0, edgeID_1, faceID_0, faceID_1;
1357c1cad2e7SMatthew G. Knepley 
13589566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(dm, e, &cone));
13599566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(bodyLabel, cone[0], &bodyID_0)); // Do I need to check the other end?
13609566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(vertexLabel, cone[0], &vertexID_0));
13619566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(vertexLabel, cone[1], &vertexID_1));
13629566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(edgeLabel, cone[0], &edgeID_0));
13639566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(edgeLabel, cone[1], &edgeID_1));
13649566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(faceLabel, cone[0], &faceID_0));
13659566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(faceLabel, cone[1], &faceID_1));
1366c1cad2e7SMatthew G. Knepley 
13679566063dSJacob Faibussowitsch     PetscCall(DMLabelSetValue(bodyLabel, e, bodyID_0));
1368c1cad2e7SMatthew G. Knepley 
13699566063dSJacob Faibussowitsch     if (edgeID_0 == edgeID_1) PetscCall(DMLabelSetValue(edgeLabel, e, edgeID_0));
13709566063dSJacob Faibussowitsch     else if (vertexID_0 > 0 && edgeID_1 > 0) PetscCall(DMLabelSetValue(edgeLabel, e, edgeID_1));
13719566063dSJacob Faibussowitsch     else if (vertexID_1 > 0 && edgeID_0 > 0) PetscCall(DMLabelSetValue(edgeLabel, e, edgeID_0));
1372c1cad2e7SMatthew G. Knepley     else { /* Do Nothing */ }
1373c1cad2e7SMatthew G. Knepley   }
1374c1cad2e7SMatthew G. Knepley 
1375c1cad2e7SMatthew G. Knepley   /* Set DMLabels for Cells */
1376c1cad2e7SMatthew G. Knepley   for (int f = fStart; f < fEnd; ++f) {
1377c1cad2e7SMatthew G. Knepley     int           edgeID_0;
1378c1cad2e7SMatthew G. Knepley     PetscInt      triBodyVal, triFaceVal;
1379c1cad2e7SMatthew G. Knepley     PetscHashIter TFLiter, TBLiter;
1380c1cad2e7SMatthew G. Knepley     PetscBool     TFLfound, TBLfound;
1381c1cad2e7SMatthew G. Knepley 
1382c1cad2e7SMatthew G. Knepley     // Convert to Hash Table
13839566063dSJacob Faibussowitsch     PetscCall(PetscHMapIFind(triFaceIDLabelMap, f - fStart, &TFLiter, &TFLfound));
138428b400f6SJacob Faibussowitsch     PetscCheck(TFLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in triFaceIDLabelMap", f);
13859566063dSJacob Faibussowitsch     PetscCall(PetscHMapIGet(triFaceIDLabelMap, f - fStart, &triFaceVal));
1386c1cad2e7SMatthew G. Knepley 
13879566063dSJacob Faibussowitsch     PetscCall(PetscHMapIFind(triBodyIDLabelMap, f - fStart, &TBLiter, &TBLfound));
138828b400f6SJacob Faibussowitsch     PetscCheck(TBLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in triBodyIDLabelMap", f);
13899566063dSJacob Faibussowitsch     PetscCall(PetscHMapIGet(triBodyIDLabelMap, f - fStart, &triBodyVal));
1390c1cad2e7SMatthew G. Knepley 
13919566063dSJacob Faibussowitsch     PetscCall(DMLabelSetValue(bodyLabel, f, triBodyVal));
13929566063dSJacob Faibussowitsch     PetscCall(DMLabelSetValue(faceLabel, f, triFaceVal));
1393c1cad2e7SMatthew G. Knepley 
1394c1cad2e7SMatthew G. Knepley     /* Finish Labeling previously unlabeled DMPlex Edges - Assumes Triangular Cell (3 Edges Max) */
13959566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(dm, f, &cone));
1396c1cad2e7SMatthew G. Knepley 
1397c1cad2e7SMatthew G. Knepley     for (int jj = 0; jj < 3; ++jj) {
13989566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValue(edgeLabel, cone[jj], &edgeID_0));
1399c1cad2e7SMatthew G. Knepley 
1400c1cad2e7SMatthew G. Knepley       if (edgeID_0 < 0) {
14019566063dSJacob Faibussowitsch         PetscCall(DMLabelSetValue(bodyLabel, cone[jj], triBodyVal));
14029566063dSJacob Faibussowitsch         PetscCall(DMLabelSetValue(faceLabel, cone[jj], triFaceVal));
1403c1cad2e7SMatthew G. Knepley       }
1404c1cad2e7SMatthew G. Knepley     }
1405c1cad2e7SMatthew G. Knepley   }
1406c1cad2e7SMatthew G. Knepley 
1407c1cad2e7SMatthew G. Knepley   *newdm = dm;
1408c1cad2e7SMatthew G. Knepley   PetscFunctionReturn(0);
1409c1cad2e7SMatthew G. Knepley }
14107bee2925SMatthew Knepley #endif
14117bee2925SMatthew Knepley 
1412c1cad2e7SMatthew G. Knepley /*@
1413c1cad2e7SMatthew 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.
1414c1cad2e7SMatthew G. Knepley 
1415c1cad2e7SMatthew G. Knepley   Collective on dm
1416c1cad2e7SMatthew G. Knepley 
1417c1cad2e7SMatthew G. Knepley   Input Parameter:
1418c1cad2e7SMatthew G. Knepley . dm - The uninflated DM object representing the mesh
1419c1cad2e7SMatthew G. Knepley 
1420c1cad2e7SMatthew G. Knepley   Output Parameter:
1421c1cad2e7SMatthew G. Knepley . dm - The inflated DM object representing the mesh
1422c1cad2e7SMatthew G. Knepley 
1423c1cad2e7SMatthew G. Knepley   Level: intermediate
1424c1cad2e7SMatthew G. Knepley 
1425db781477SPatrick Sanan .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateEGADS()`
1426c1cad2e7SMatthew G. Knepley @*/
1427*9371c9d4SSatish Balay PetscErrorCode DMPlexInflateToGeomModel(DM dm) {
1428c1cad2e7SMatthew G. Knepley #if defined(PETSC_HAVE_EGADS)
1429c1cad2e7SMatthew G. Knepley   /* EGADS Variables */
1430c1cad2e7SMatthew G. Knepley   ego            model, geom, body, face, edge;
1431c1cad2e7SMatthew G. Knepley   ego           *bodies;
1432c1cad2e7SMatthew G. Knepley   int            Nb, oclass, mtype, *senses;
1433c1cad2e7SMatthew G. Knepley   double         result[3];
1434c1cad2e7SMatthew G. Knepley   /* PETSc Variables */
1435c1cad2e7SMatthew G. Knepley   DM             cdm;
1436c1cad2e7SMatthew G. Knepley   PetscContainer modelObj;
1437c1cad2e7SMatthew G. Knepley   DMLabel        bodyLabel, faceLabel, edgeLabel, vertexLabel;
1438c1cad2e7SMatthew G. Knepley   Vec            coordinates;
1439c1cad2e7SMatthew G. Knepley   PetscScalar   *coords;
1440c1cad2e7SMatthew G. Knepley   PetscInt       bodyID, faceID, edgeID, vertexID;
1441c1cad2e7SMatthew G. Knepley   PetscInt       cdim, d, vStart, vEnd, v;
1442c1cad2e7SMatthew G. Knepley #endif
1443c1cad2e7SMatthew G. Knepley 
1444c1cad2e7SMatthew G. Knepley   PetscFunctionBegin;
1445c1cad2e7SMatthew G. Knepley #if defined(PETSC_HAVE_EGADS)
14469566063dSJacob Faibussowitsch   PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj));
1447c1cad2e7SMatthew G. Knepley   if (!modelObj) PetscFunctionReturn(0);
14489566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
14499566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
14509566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
14519566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel));
14529566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel));
14539566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel));
14549566063dSJacob Faibussowitsch   PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel));
1455c1cad2e7SMatthew G. Knepley 
14569566063dSJacob Faibussowitsch   PetscCall(PetscContainerGetPointer(modelObj, (void **)&model));
14579566063dSJacob Faibussowitsch   PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses));
1458c1cad2e7SMatthew G. Knepley 
14599566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
14609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(coordinates, &coords));
1461c1cad2e7SMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
1462c1cad2e7SMatthew G. Knepley     PetscScalar *vcoords;
1463c1cad2e7SMatthew G. Knepley 
14649566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(bodyLabel, v, &bodyID));
14659566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(faceLabel, v, &faceID));
14669566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(edgeLabel, v, &edgeID));
14679566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(vertexLabel, v, &vertexID));
1468c1cad2e7SMatthew G. Knepley 
146963a3b9bcSJacob Faibussowitsch     PetscCheck(bodyID < Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %" PetscInt_FMT " is not in [0, %d)", bodyID, Nb);
1470c1cad2e7SMatthew G. Knepley     body = bodies[bodyID];
1471c1cad2e7SMatthew G. Knepley 
14729566063dSJacob Faibussowitsch     PetscCall(DMPlexPointLocalRef(cdm, v, coords, (void *)&vcoords));
1473c1cad2e7SMatthew G. Knepley     if (edgeID > 0) {
1474c1cad2e7SMatthew G. Knepley       /* Snap to EDGE at nearest location */
1475c1cad2e7SMatthew G. Knepley       double params[1];
14769566063dSJacob Faibussowitsch       PetscCall(EG_objectBodyTopo(body, EDGE, edgeID, &edge));
14779566063dSJacob Faibussowitsch       PetscCall(EG_invEvaluate(edge, vcoords, params, result)); // Get (x,y,z) of nearest point on EDGE
1478c1cad2e7SMatthew G. Knepley       for (d = 0; d < cdim; ++d) vcoords[d] = result[d];
1479c1cad2e7SMatthew G. Knepley     } else if (faceID > 0) {
1480c1cad2e7SMatthew G. Knepley       /* Snap to FACE at nearest location */
1481c1cad2e7SMatthew G. Knepley       double params[2];
14829566063dSJacob Faibussowitsch       PetscCall(EG_objectBodyTopo(body, FACE, faceID, &face));
14839566063dSJacob Faibussowitsch       PetscCall(EG_invEvaluate(face, vcoords, params, result)); // Get (x,y,z) of nearest point on FACE
1484c1cad2e7SMatthew G. Knepley       for (d = 0; d < cdim; ++d) vcoords[d] = result[d];
1485c1cad2e7SMatthew G. Knepley     }
1486c1cad2e7SMatthew G. Knepley   }
14879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(coordinates, &coords));
1488c1cad2e7SMatthew G. Knepley   /* Clear out global coordinates */
14896858538eSMatthew G. Knepley   PetscCall(VecDestroy(&dm->coordinates[0].x));
1490c1cad2e7SMatthew G. Knepley #endif
1491c1cad2e7SMatthew G. Knepley   PetscFunctionReturn(0);
1492c1cad2e7SMatthew G. Knepley }
1493c1cad2e7SMatthew G. Knepley 
14947bee2925SMatthew Knepley /*@C
1495c1cad2e7SMatthew G. Knepley   DMPlexCreateEGADSFromFile - Create a DMPlex mesh from an EGADS, IGES, or STEP file.
14967bee2925SMatthew Knepley 
14977bee2925SMatthew Knepley   Collective
14987bee2925SMatthew Knepley 
14997bee2925SMatthew Knepley   Input Parameters:
15007bee2925SMatthew Knepley + comm     - The MPI communicator
1501c1cad2e7SMatthew G. Knepley - filename - The name of the EGADS, IGES, or STEP file
15027bee2925SMatthew Knepley 
15037bee2925SMatthew Knepley   Output Parameter:
15047bee2925SMatthew Knepley . dm       - The DM object representing the mesh
15057bee2925SMatthew Knepley 
15067bee2925SMatthew Knepley   Level: beginner
15077bee2925SMatthew Knepley 
1508db781477SPatrick Sanan .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateEGADS()`, `DMPlexCreateEGADSLiteFromFile()`
15097bee2925SMatthew Knepley @*/
1510*9371c9d4SSatish Balay PetscErrorCode DMPlexCreateEGADSFromFile(MPI_Comm comm, const char filename[], DM *dm) {
15117bee2925SMatthew Knepley   PetscMPIInt rank;
15127bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS)
15137bee2925SMatthew Knepley   ego context = NULL, model = NULL;
15147bee2925SMatthew Knepley #endif
1515c1cad2e7SMatthew G. Knepley   PetscBool printModel = PETSC_FALSE, tessModel = PETSC_FALSE, newModel = PETSC_FALSE;
15167bee2925SMatthew Knepley 
15177bee2925SMatthew Knepley   PetscFunctionBegin;
15187bee2925SMatthew Knepley   PetscValidCharPointer(filename, 2);
15199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_print_model", &printModel, NULL));
15209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_tess_model", &tessModel, NULL));
15219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_new_model", &newModel, NULL));
15229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
15237bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS)
1524dd400576SPatrick Sanan   if (rank == 0) {
15259566063dSJacob Faibussowitsch     PetscCall(EG_open(&context));
15269566063dSJacob Faibussowitsch     PetscCall(EG_loadModel(context, 0, filename, &model));
15279566063dSJacob Faibussowitsch     if (printModel) PetscCall(DMPlexEGADSPrintModel_Internal(model));
15287bee2925SMatthew Knepley   }
15299566063dSJacob Faibussowitsch   if (tessModel) PetscCall(DMPlexCreateEGADS_Tess_Internal(comm, context, model, dm));
15309566063dSJacob Faibussowitsch   else if (newModel) PetscCall(DMPlexCreateEGADS(comm, context, model, dm));
15319566063dSJacob Faibussowitsch   else PetscCall(DMPlexCreateEGADS_Internal(comm, context, model, dm));
1532c03d1177SJose E. Roman   PetscFunctionReturn(0);
15337bee2925SMatthew Knepley #else
1534c1cad2e7SMatthew G. Knepley   SETERRQ(comm, PETSC_ERR_SUP, "This method requires EGADS support. Reconfigure using --download-egads");
15357bee2925SMatthew Knepley #endif
15367bee2925SMatthew Knepley }
1537