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 25d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSnapToGeomModel_EGADS_Internal(DM dm, PetscInt p, ego model, PetscInt bodyID, PetscInt faceID, PetscInt edgeID, const PetscScalar mcoords[], PetscScalar gcoords[]) 26d71ae5a4SJacob Faibussowitsch { 27c1cad2e7SMatthew G. Knepley DM cdm; 28c1cad2e7SMatthew G. Knepley ego *bodies; 29c1cad2e7SMatthew G. Knepley ego geom, body, obj; 30a5b23f4aSJose E. Roman /* result has to hold derivatives, along with the value */ 31c1cad2e7SMatthew G. Knepley double params[3], result[18], paramsV[16 * 3], resultV[16 * 3], range[4]; 32c1cad2e7SMatthew G. Knepley int Nb, oclass, mtype, *senses, peri; 33c1cad2e7SMatthew G. Knepley Vec coordinatesLocal; 34c1cad2e7SMatthew G. Knepley PetscScalar *coords = NULL; 35c1cad2e7SMatthew G. Knepley PetscInt Nv, v, Np = 0, pm; 36c1cad2e7SMatthew G. Knepley PetscInt dE, d; 37c1cad2e7SMatthew G. Knepley 38c1cad2e7SMatthew G. Knepley PetscFunctionBeginHot; 399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dE)); 419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinatesLocal)); 429566063dSJacob Faibussowitsch PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 4363a3b9bcSJacob Faibussowitsch PetscCheck(bodyID < Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %" PetscInt_FMT " is not in [0, %d)", bodyID, Nb); 44c1cad2e7SMatthew G. Knepley body = bodies[bodyID]; 45c1cad2e7SMatthew G. Knepley 469371c9d4SSatish Balay if (edgeID >= 0) { 479371c9d4SSatish Balay PetscCall(EG_objectBodyTopo(body, EDGE, edgeID, &obj)); 489371c9d4SSatish Balay Np = 1; 499371c9d4SSatish Balay } else if (faceID >= 0) { 509371c9d4SSatish Balay PetscCall(EG_objectBodyTopo(body, FACE, faceID, &obj)); 519371c9d4SSatish Balay Np = 2; 529371c9d4SSatish Balay } else { 53c1cad2e7SMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 54c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 55c1cad2e7SMatthew G. Knepley } 56c1cad2e7SMatthew G. Knepley 57c1cad2e7SMatthew G. Knepley /* Calculate parameters (t or u,v) for vertices */ 589566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords)); 59c1cad2e7SMatthew G. Knepley Nv /= dE; 60c1cad2e7SMatthew G. Knepley if (Nv == 1) { 619566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords)); 62c1cad2e7SMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 63c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 64c1cad2e7SMatthew G. Knepley } 6563a3b9bcSJacob Faibussowitsch PetscCheck(Nv <= 16, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %" PetscInt_FMT " coordinates associated to point %" PetscInt_FMT, Nv, p); 66c1cad2e7SMatthew G. Knepley 67c1cad2e7SMatthew G. Knepley /* Correct EGADSlite 2pi bug when calculating nearest point on Periodic Surfaces */ 689566063dSJacob Faibussowitsch PetscCall(EG_getRange(obj, range, &peri)); 69c1cad2e7SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 709566063dSJacob Faibussowitsch PetscCall(EG_invEvaluate(obj, &coords[v * dE], ¶msV[v * 3], &resultV[v * 3])); 71c1cad2e7SMatthew G. Knepley #if 1 72c1cad2e7SMatthew G. Knepley if (peri > 0) { 739371c9d4SSatish Balay if (paramsV[v * 3 + 0] + 1.e-4 < range[0]) { 749371c9d4SSatish Balay paramsV[v * 3 + 0] += 2. * PETSC_PI; 759371c9d4SSatish Balay } else if (paramsV[v * 3 + 0] - 1.e-4 > range[1]) { 769371c9d4SSatish Balay paramsV[v * 3 + 0] -= 2. * PETSC_PI; 779371c9d4SSatish Balay } 78c1cad2e7SMatthew G. Knepley } 79c1cad2e7SMatthew G. Knepley if (peri > 1) { 809371c9d4SSatish Balay if (paramsV[v * 3 + 1] + 1.e-4 < range[2]) { 819371c9d4SSatish Balay paramsV[v * 3 + 1] += 2. * PETSC_PI; 829371c9d4SSatish Balay } else if (paramsV[v * 3 + 1] - 1.e-4 > range[3]) { 839371c9d4SSatish Balay paramsV[v * 3 + 1] -= 2. * PETSC_PI; 849371c9d4SSatish Balay } 85c1cad2e7SMatthew G. Knepley } 86c1cad2e7SMatthew G. Knepley #endif 87c1cad2e7SMatthew G. Knepley } 889566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords)); 89c1cad2e7SMatthew G. Knepley /* Calculate parameters (t or u,v) for new vertex at edge midpoint */ 90c1cad2e7SMatthew G. Knepley for (pm = 0; pm < Np; ++pm) { 91c1cad2e7SMatthew G. Knepley params[pm] = 0.; 92c1cad2e7SMatthew G. Knepley for (v = 0; v < Nv; ++v) params[pm] += paramsV[v * 3 + pm]; 93c1cad2e7SMatthew G. Knepley params[pm] /= Nv; 94c1cad2e7SMatthew G. Knepley } 9563a3b9bcSJacob Faibussowitsch PetscCheck(!(params[0] < range[0]) && !(params[0] > range[1]), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " had bad interpolation", p); 961dca8a05SBarry 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); 97c1cad2e7SMatthew G. Knepley /* Put coordinates for new vertex in result[] */ 989566063dSJacob Faibussowitsch PetscCall(EG_evaluate(obj, params, result)); 99c1cad2e7SMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = result[d]; 100c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 101c1cad2e7SMatthew G. Knepley } 102c1cad2e7SMatthew G. Knepley #endif 103c1cad2e7SMatthew G. Knepley 104a8ededdfSMatthew G. Knepley /*@ 105a8ededdfSMatthew 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. 106a8ededdfSMatthew G. Knepley 107a8ededdfSMatthew G. Knepley Not collective 108a8ededdfSMatthew G. Knepley 109a8ededdfSMatthew G. Knepley Input Parameters: 110a8ededdfSMatthew G. Knepley + dm - The DMPlex object 111a8ededdfSMatthew G. Knepley . p - The mesh point 112d410b0cfSMatthew G. Knepley . dE - The coordinate dimension 113a8ededdfSMatthew G. Knepley - mcoords - A coordinate point lying on the mesh point 114a8ededdfSMatthew G. Knepley 115a8ededdfSMatthew G. Knepley Output Parameter: 116a8ededdfSMatthew G. Knepley . gcoords - The closest coordinate point on the geometry model associated with 'p' to the given point 117a8ededdfSMatthew G. Knepley 118d410b0cfSMatthew 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. 119a8ededdfSMatthew G. Knepley 120a8ededdfSMatthew G. Knepley Level: intermediate 121a8ededdfSMatthew G. Knepley 122db781477SPatrick Sanan .seealso: `DMRefine()`, `DMPlexCreate()`, `DMPlexSetRefinementUniform()` 123a8ededdfSMatthew G. Knepley @*/ 124d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexSnapToGeomModel(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[]) 125d71ae5a4SJacob Faibussowitsch { 126d410b0cfSMatthew G. Knepley PetscInt d; 127a8ededdfSMatthew G. Knepley 128c1cad2e7SMatthew G. Knepley PetscFunctionBeginHot; 129a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS 130c1cad2e7SMatthew G. Knepley { 131c1cad2e7SMatthew G. Knepley DM_Plex *plex = (DM_Plex *)dm->data; 132c1cad2e7SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel; 133c1cad2e7SMatthew G. Knepley PetscInt bodyID, faceID, edgeID; 134c1cad2e7SMatthew G. Knepley PetscContainer modelObj; 135c1cad2e7SMatthew G. Knepley ego model; 136c1cad2e7SMatthew G. Knepley PetscBool islite = PETSC_FALSE; 137c1cad2e7SMatthew G. Knepley 1389566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 1399566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 1409566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel)); 141c1cad2e7SMatthew G. Knepley if (!bodyLabel || !faceLabel || !edgeLabel || plex->ignoreModel) { 142f0fcf0adSMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 143a8ededdfSMatthew G. Knepley PetscFunctionReturn(0); 144a8ededdfSMatthew G. Knepley } 1459566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj)); 146c1cad2e7SMatthew G. Knepley if (!modelObj) { 1479566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "EGADSLite Model", (PetscObject *)&modelObj)); 148c1cad2e7SMatthew G. Knepley islite = PETSC_TRUE; 149c1cad2e7SMatthew G. Knepley } 150c1cad2e7SMatthew G. Knepley if (!modelObj) { 151c1cad2e7SMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 152c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 153c1cad2e7SMatthew G. Knepley } 1549566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(modelObj, (void **)&model)); 1559566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, p, &bodyID)); 1569566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(faceLabel, p, &faceID)); 1579566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(edgeLabel, p, &edgeID)); 158c1cad2e7SMatthew G. Knepley /* Allows for "Connective" Plex Edges present in models with multiple non-touching Entities */ 159c1cad2e7SMatthew G. Knepley if (bodyID < 0) { 160f0fcf0adSMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 161f0fcf0adSMatthew G. Knepley PetscFunctionReturn(0); 162a8ededdfSMatthew G. Knepley } 1639566063dSJacob Faibussowitsch if (islite) PetscCall(DMPlexSnapToGeomModel_EGADSLite_Internal(dm, p, model, bodyID, faceID, edgeID, mcoords, gcoords)); 1649566063dSJacob Faibussowitsch else PetscCall(DMPlexSnapToGeomModel_EGADS_Internal(dm, p, model, bodyID, faceID, edgeID, mcoords, gcoords)); 165f0fcf0adSMatthew G. Knepley } 166a8ededdfSMatthew G. Knepley #else 167f0fcf0adSMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 168a8ededdfSMatthew G. Knepley #endif 169a8ededdfSMatthew G. Knepley PetscFunctionReturn(0); 170a8ededdfSMatthew G. Knepley } 1717bee2925SMatthew Knepley 1727bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS) 173d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexEGADSPrintModel_Internal(ego model) 174d71ae5a4SJacob Faibussowitsch { 1757bee2925SMatthew Knepley ego geom, *bodies, *objs, *nobjs, *mobjs, *lobjs; 1767bee2925SMatthew Knepley int oclass, mtype, *senses; 1777bee2925SMatthew Knepley int Nb, b; 1787bee2925SMatthew Knepley 1797bee2925SMatthew Knepley PetscFunctionBeginUser; 1807bee2925SMatthew Knepley /* test bodyTopo functions */ 1819566063dSJacob Faibussowitsch PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 1829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of BODIES (nbodies): %d \n", Nb)); 1837bee2925SMatthew Knepley 1847bee2925SMatthew Knepley for (b = 0; b < Nb; ++b) { 1857bee2925SMatthew Knepley ego body = bodies[b]; 1867bee2925SMatthew Knepley int id, Nsh, Nf, Nl, l, Ne, e, Nv, v; 1877bee2925SMatthew Knepley 1887bee2925SMatthew Knepley /* Output Basic Model Topology */ 1899566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, SHELL, &Nsh, &objs)); 1909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of SHELLS: %d \n", Nsh)); 1917bee2925SMatthew Knepley EG_free(objs); 1927bee2925SMatthew Knepley 1939566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &objs)); 1949566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of FACES: %d \n", Nf)); 1957bee2925SMatthew Knepley EG_free(objs); 1967bee2925SMatthew Knepley 1979566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs)); 1989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of LOOPS: %d \n", Nl)); 1997bee2925SMatthew Knepley 2009566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &objs)); 2019566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of EDGES: %d \n", Ne)); 2027bee2925SMatthew Knepley EG_free(objs); 2037bee2925SMatthew Knepley 2049566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &objs)); 2059566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Number of NODES: %d \n", Nv)); 2067bee2925SMatthew Knepley EG_free(objs); 2077bee2925SMatthew Knepley 2087bee2925SMatthew Knepley for (l = 0; l < Nl; ++l) { 2097bee2925SMatthew Knepley ego loop = lobjs[l]; 2107bee2925SMatthew Knepley 2117bee2925SMatthew Knepley id = EG_indexBodyTopo(body, loop); 2129566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " LOOP ID: %d\n", id)); 2137bee2925SMatthew Knepley 2147bee2925SMatthew Knepley /* Get EDGE info which associated with the current LOOP */ 2159566063dSJacob Faibussowitsch PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses)); 2167bee2925SMatthew Knepley 2177bee2925SMatthew Knepley for (e = 0; e < Ne; ++e) { 2187bee2925SMatthew Knepley ego edge = objs[e]; 2197bee2925SMatthew Knepley double range[4] = {0., 0., 0., 0.}; 2207bee2925SMatthew Knepley double point[3] = {0., 0., 0.}; 221266cfabeSMatthew G. Knepley double params[3] = {0., 0., 0.}; 222266cfabeSMatthew G. Knepley double result[18]; 2237bee2925SMatthew Knepley int peri; 2247bee2925SMatthew Knepley 2255f80ce2aSJacob Faibussowitsch id = EG_indexBodyTopo(body, edge); 2269566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " EDGE ID: %d (%d)\n", id, e)); 2277bee2925SMatthew Knepley 2289566063dSJacob Faibussowitsch PetscCall(EG_getRange(edge, range, &peri)); 2299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " Range = %lf, %lf, %lf, %lf \n", range[0], range[1], range[2], range[3])); 2307bee2925SMatthew Knepley 2317bee2925SMatthew Knepley /* Get NODE info which associated with the current EDGE */ 2329566063dSJacob Faibussowitsch PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 233266cfabeSMatthew G. Knepley if (mtype == DEGENERATE) { 2349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " EDGE %d is DEGENERATE \n", id)); 235266cfabeSMatthew G. Knepley } else { 236266cfabeSMatthew G. Knepley params[0] = range[0]; 2379566063dSJacob Faibussowitsch PetscCall(EG_evaluate(edge, params, result)); 2389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " between (%lf, %lf, %lf)", result[0], result[1], result[2])); 239266cfabeSMatthew G. Knepley params[0] = range[1]; 2409566063dSJacob Faibussowitsch PetscCall(EG_evaluate(edge, params, result)); 2419566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " and (%lf, %lf, %lf)\n", result[0], result[1], result[2])); 242266cfabeSMatthew G. Knepley } 2437bee2925SMatthew Knepley 2447bee2925SMatthew Knepley for (v = 0; v < Nv; ++v) { 2457bee2925SMatthew Knepley ego vertex = nobjs[v]; 2467bee2925SMatthew Knepley double limits[4]; 2477bee2925SMatthew Knepley int dummy; 2487bee2925SMatthew Knepley 2499566063dSJacob Faibussowitsch PetscCall(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses)); 2507bee2925SMatthew Knepley id = EG_indexBodyTopo(body, vertex); 2519566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " NODE ID: %d \n", id)); 2529566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " (x, y, z) = (%lf, %lf, %lf) \n", limits[0], limits[1], limits[2])); 2537bee2925SMatthew Knepley 2547bee2925SMatthew Knepley point[0] = point[0] + limits[0]; 2557bee2925SMatthew Knepley point[1] = point[1] + limits[1]; 2567bee2925SMatthew Knepley point[2] = point[2] + limits[2]; 2577bee2925SMatthew Knepley } 2587bee2925SMatthew Knepley } 2597bee2925SMatthew Knepley } 2607bee2925SMatthew Knepley EG_free(lobjs); 2617bee2925SMatthew Knepley } 2627bee2925SMatthew Knepley PetscFunctionReturn(0); 2637bee2925SMatthew Knepley } 2647bee2925SMatthew Knepley 265d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexEGADSDestroy_Private(void *context) 266d71ae5a4SJacob Faibussowitsch { 2677bee2925SMatthew Knepley if (context) EG_close((ego)context); 2687bee2925SMatthew Knepley return 0; 2697bee2925SMatthew Knepley } 2707bee2925SMatthew Knepley 271d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateEGADS_Internal(MPI_Comm comm, ego context, ego model, DM *newdm) 272d71ae5a4SJacob Faibussowitsch { 273c1cad2e7SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel; 2747bee2925SMatthew Knepley PetscInt cStart, cEnd, c; 2757bee2925SMatthew Knepley /* EGADSLite variables */ 2767bee2925SMatthew Knepley ego geom, *bodies, *objs, *nobjs, *mobjs, *lobjs; 2777bee2925SMatthew Knepley int oclass, mtype, nbodies, *senses; 2787bee2925SMatthew Knepley int b; 2797bee2925SMatthew Knepley /* PETSc variables */ 2807bee2925SMatthew Knepley DM dm; 2817bee2925SMatthew Knepley PetscHMapI edgeMap = NULL; 2827bee2925SMatthew 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; 2837bee2925SMatthew Knepley PetscInt *cells = NULL, *cone = NULL; 2847bee2925SMatthew Knepley PetscReal *coords = NULL; 2857bee2925SMatthew Knepley PetscMPIInt rank; 2867bee2925SMatthew Knepley 2877bee2925SMatthew Knepley PetscFunctionBegin; 2889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 289dd400576SPatrick Sanan if (rank == 0) { 290266cfabeSMatthew G. Knepley const PetscInt debug = 0; 291266cfabeSMatthew G. Knepley 2927bee2925SMatthew Knepley /* --------------------------------------------------------------------------------------------------- 2937bee2925SMatthew Knepley Generate Petsc Plex 2947bee2925SMatthew Knepley Get all Nodes in model, record coordinates in a correctly formatted array 2957bee2925SMatthew Knepley Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array 2967bee2925SMatthew Knepley We need to uniformly refine the initial geometry to guarantee a valid mesh 2977bee2925SMatthew Knepley */ 2987bee2925SMatthew Knepley 2997bee2925SMatthew Knepley /* Calculate cell and vertex sizes */ 3009566063dSJacob Faibussowitsch PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses)); 3019566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&edgeMap)); 3027bee2925SMatthew Knepley numEdges = 0; 3037bee2925SMatthew Knepley for (b = 0; b < nbodies; ++b) { 3047bee2925SMatthew Knepley ego body = bodies[b]; 3057bee2925SMatthew Knepley int id, Nl, l, Nv, v; 3067bee2925SMatthew Knepley 3079566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs)); 3087bee2925SMatthew Knepley for (l = 0; l < Nl; ++l) { 3097bee2925SMatthew Knepley ego loop = lobjs[l]; 3107bee2925SMatthew Knepley int Ner = 0, Ne, e, Nc; 3117bee2925SMatthew Knepley 3129566063dSJacob Faibussowitsch PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses)); 3137bee2925SMatthew Knepley for (e = 0; e < Ne; ++e) { 3147bee2925SMatthew Knepley ego edge = objs[e]; 3157bee2925SMatthew Knepley int Nv, id; 3167bee2925SMatthew Knepley PetscHashIter iter; 3177bee2925SMatthew Knepley PetscBool found; 3187bee2925SMatthew Knepley 3199566063dSJacob Faibussowitsch PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 3207bee2925SMatthew Knepley if (mtype == DEGENERATE) continue; 3215f80ce2aSJacob Faibussowitsch id = EG_indexBodyTopo(body, edge); 3229566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(edgeMap, id - 1, &iter, &found)); 3239566063dSJacob Faibussowitsch if (!found) PetscCall(PetscHMapISet(edgeMap, id - 1, numEdges++)); 3247bee2925SMatthew Knepley ++Ner; 3257bee2925SMatthew Knepley } 3269371c9d4SSatish Balay if (Ner == 2) { 3279371c9d4SSatish Balay Nc = 2; 3289371c9d4SSatish Balay } else if (Ner == 3) { 3299371c9d4SSatish Balay Nc = 4; 3309371c9d4SSatish Balay } else if (Ner == 4) { 3319371c9d4SSatish Balay Nc = 8; 3329371c9d4SSatish Balay ++numQuads; 3339371c9d4SSatish Balay } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot support loop with %d edges", Ner); 3347bee2925SMatthew Knepley numCells += Nc; 3357bee2925SMatthew Knepley newCells += Nc - 1; 3367bee2925SMatthew Knepley maxCorners = PetscMax(Ner * 2 + 1, maxCorners); 3377bee2925SMatthew Knepley } 3389566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs)); 3397bee2925SMatthew Knepley for (v = 0; v < Nv; ++v) { 3407bee2925SMatthew Knepley ego vertex = nobjs[v]; 3417bee2925SMatthew Knepley 3427bee2925SMatthew Knepley id = EG_indexBodyTopo(body, vertex); 3437bee2925SMatthew Knepley /* TODO: Instead of assuming contiguous ids, we could use a hash table */ 3447bee2925SMatthew Knepley numVertices = PetscMax(id, numVertices); 3457bee2925SMatthew Knepley } 3467bee2925SMatthew Knepley EG_free(lobjs); 3477bee2925SMatthew Knepley EG_free(nobjs); 3487bee2925SMatthew Knepley } 3499566063dSJacob Faibussowitsch PetscCall(PetscHMapIGetSize(edgeMap, &numEdges)); 3507bee2925SMatthew Knepley newVertices = numEdges + numQuads; 3517bee2925SMatthew Knepley numVertices += newVertices; 3527bee2925SMatthew Knepley 3537bee2925SMatthew Knepley dim = 2; /* Assume 3D Models :: Need to update to handle 2D Models in the future */ 3547bee2925SMatthew Knepley cdim = 3; /* Assume 3D Models :: Need to update to handle 2D Models in the future */ 3557bee2925SMatthew Knepley numCorners = 3; /* Split cells into triangles */ 3569566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(numVertices * cdim, &coords, numCells * numCorners, &cells, maxCorners, &cone)); 3577bee2925SMatthew Knepley 3587bee2925SMatthew Knepley /* Get vertex coordinates */ 3597bee2925SMatthew Knepley for (b = 0; b < nbodies; ++b) { 3607bee2925SMatthew Knepley ego body = bodies[b]; 3617bee2925SMatthew Knepley int id, Nv, v; 3627bee2925SMatthew Knepley 3639566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs)); 3647bee2925SMatthew Knepley for (v = 0; v < Nv; ++v) { 3657bee2925SMatthew Knepley ego vertex = nobjs[v]; 3667bee2925SMatthew Knepley double limits[4]; 3677bee2925SMatthew Knepley int dummy; 3687bee2925SMatthew Knepley 3699566063dSJacob Faibussowitsch PetscCall(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses)); 3705f80ce2aSJacob Faibussowitsch id = EG_indexBodyTopo(body, vertex); 3717bee2925SMatthew Knepley coords[(id - 1) * cdim + 0] = limits[0]; 3727bee2925SMatthew Knepley coords[(id - 1) * cdim + 1] = limits[1]; 3737bee2925SMatthew Knepley coords[(id - 1) * cdim + 2] = limits[2]; 3747bee2925SMatthew Knepley } 3757bee2925SMatthew Knepley EG_free(nobjs); 3767bee2925SMatthew Knepley } 3779566063dSJacob Faibussowitsch PetscCall(PetscHMapIClear(edgeMap)); 3787bee2925SMatthew Knepley fOff = numVertices - newVertices + numEdges; 3797bee2925SMatthew Knepley numEdges = 0; 3807bee2925SMatthew Knepley numQuads = 0; 3817bee2925SMatthew Knepley for (b = 0; b < nbodies; ++b) { 3827bee2925SMatthew Knepley ego body = bodies[b]; 3837bee2925SMatthew Knepley int Nl, l; 3847bee2925SMatthew Knepley 3859566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs)); 3867bee2925SMatthew Knepley for (l = 0; l < Nl; ++l) { 3877bee2925SMatthew Knepley ego loop = lobjs[l]; 3887bee2925SMatthew Knepley int lid, Ner = 0, Ne, e; 3897bee2925SMatthew Knepley 3905f80ce2aSJacob Faibussowitsch lid = EG_indexBodyTopo(body, loop); 3919566063dSJacob Faibussowitsch PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses)); 3927bee2925SMatthew Knepley for (e = 0; e < Ne; ++e) { 3937bee2925SMatthew Knepley ego edge = objs[e]; 3947bee2925SMatthew Knepley int eid, Nv; 3957bee2925SMatthew Knepley PetscHashIter iter; 3967bee2925SMatthew Knepley PetscBool found; 3977bee2925SMatthew Knepley 3989566063dSJacob Faibussowitsch PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 3997bee2925SMatthew Knepley if (mtype == DEGENERATE) continue; 4007bee2925SMatthew Knepley ++Ner; 4015f80ce2aSJacob Faibussowitsch eid = EG_indexBodyTopo(body, edge); 4029566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(edgeMap, eid - 1, &iter, &found)); 4037bee2925SMatthew Knepley if (!found) { 4047bee2925SMatthew Knepley PetscInt v = numVertices - newVertices + numEdges; 4057bee2925SMatthew Knepley double range[4], params[3] = {0., 0., 0.}, result[18]; 4067bee2925SMatthew Knepley int periodic[2]; 4077bee2925SMatthew Knepley 4089566063dSJacob Faibussowitsch PetscCall(PetscHMapISet(edgeMap, eid - 1, numEdges++)); 4099566063dSJacob Faibussowitsch PetscCall(EG_getRange(edge, range, periodic)); 4107bee2925SMatthew Knepley params[0] = 0.5 * (range[0] + range[1]); 4119566063dSJacob Faibussowitsch PetscCall(EG_evaluate(edge, params, result)); 4127bee2925SMatthew Knepley coords[v * cdim + 0] = result[0]; 4137bee2925SMatthew Knepley coords[v * cdim + 1] = result[1]; 4147bee2925SMatthew Knepley coords[v * cdim + 2] = result[2]; 4157bee2925SMatthew Knepley } 4167bee2925SMatthew Knepley } 4177bee2925SMatthew Knepley if (Ner == 4) { 4187bee2925SMatthew Knepley PetscInt v = fOff + numQuads++; 419266cfabeSMatthew G. Knepley ego *fobjs, face; 4207bee2925SMatthew Knepley double range[4], params[3] = {0., 0., 0.}, result[18]; 421266cfabeSMatthew G. Knepley int Nf, fid, periodic[2]; 4227bee2925SMatthew Knepley 4239566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs)); 424266cfabeSMatthew G. Knepley face = fobjs[0]; 4255f80ce2aSJacob Faibussowitsch fid = EG_indexBodyTopo(body, face); 42608401ef6SPierre Jolivet PetscCheck(Nf == 1, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Loop %d has %d faces, instead of 1 (%d)", lid - 1, Nf, fid); 4279566063dSJacob Faibussowitsch PetscCall(EG_getRange(face, range, periodic)); 4287bee2925SMatthew Knepley params[0] = 0.5 * (range[0] + range[1]); 4297bee2925SMatthew Knepley params[1] = 0.5 * (range[2] + range[3]); 4309566063dSJacob Faibussowitsch PetscCall(EG_evaluate(face, params, result)); 4317bee2925SMatthew Knepley coords[v * cdim + 0] = result[0]; 4327bee2925SMatthew Knepley coords[v * cdim + 1] = result[1]; 4337bee2925SMatthew Knepley coords[v * cdim + 2] = result[2]; 4347bee2925SMatthew Knepley } 4357bee2925SMatthew Knepley } 4367bee2925SMatthew Knepley } 4371dca8a05SBarry Smith PetscCheck(numEdges + numQuads == newVertices, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of new vertices %" PetscInt_FMT " != %" PetscInt_FMT " previous count", numEdges + numQuads, newVertices); 4387bee2925SMatthew Knepley 4397bee2925SMatthew Knepley /* Get cell vertices by traversing loops */ 4407bee2925SMatthew Knepley numQuads = 0; 4417bee2925SMatthew Knepley cOff = 0; 4427bee2925SMatthew Knepley for (b = 0; b < nbodies; ++b) { 4437bee2925SMatthew Knepley ego body = bodies[b]; 4447bee2925SMatthew Knepley int id, Nl, l; 4457bee2925SMatthew Knepley 4469566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs)); 4477bee2925SMatthew Knepley for (l = 0; l < Nl; ++l) { 4487bee2925SMatthew Knepley ego loop = lobjs[l]; 4497bee2925SMatthew Knepley int lid, Ner = 0, Ne, e, nc = 0, c, Nt, t; 4507bee2925SMatthew Knepley 4515f80ce2aSJacob Faibussowitsch lid = EG_indexBodyTopo(body, loop); 4529566063dSJacob Faibussowitsch PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses)); 4537bee2925SMatthew Knepley 4547bee2925SMatthew Knepley for (e = 0; e < Ne; ++e) { 4557bee2925SMatthew Knepley ego edge = objs[e]; 4567bee2925SMatthew Knepley int points[3]; 4577bee2925SMatthew Knepley int eid, Nv, v, tmp; 4587bee2925SMatthew Knepley 4597bee2925SMatthew Knepley eid = EG_indexBodyTopo(body, edge); 4609566063dSJacob Faibussowitsch PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 461266cfabeSMatthew G. Knepley if (mtype == DEGENERATE) continue; 462266cfabeSMatthew G. Knepley else ++Ner; 46308401ef6SPierre Jolivet PetscCheck(Nv == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Edge %d has %d vertices != 2", eid, Nv); 4647bee2925SMatthew Knepley 4657bee2925SMatthew Knepley for (v = 0; v < Nv; ++v) { 4667bee2925SMatthew Knepley ego vertex = nobjs[v]; 4677bee2925SMatthew Knepley 4687bee2925SMatthew Knepley id = EG_indexBodyTopo(body, vertex); 4697bee2925SMatthew Knepley points[v * 2] = id - 1; 4707bee2925SMatthew Knepley } 4717bee2925SMatthew Knepley { 4727bee2925SMatthew Knepley PetscInt edgeNum; 4737bee2925SMatthew Knepley 4749566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(edgeMap, eid - 1, &edgeNum)); 4757bee2925SMatthew Knepley points[1] = numVertices - newVertices + edgeNum; 4767bee2925SMatthew Knepley } 4777bee2925SMatthew Knepley /* EGADS loops are not oriented, but seem to be in order, so we must piece them together */ 4787bee2925SMatthew Knepley if (!nc) { 4797bee2925SMatthew Knepley for (v = 0; v < Nv + 1; ++v) cone[nc++] = points[v]; 4807bee2925SMatthew Knepley } else { 4819371c9d4SSatish Balay if (cone[nc - 1] == points[0]) { 4829371c9d4SSatish Balay cone[nc++] = points[1]; 4839371c9d4SSatish Balay if (cone[0] != points[2]) cone[nc++] = points[2]; 4849371c9d4SSatish Balay } else if (cone[nc - 1] == points[2]) { 4859371c9d4SSatish Balay cone[nc++] = points[1]; 4869371c9d4SSatish Balay if (cone[0] != points[0]) cone[nc++] = points[0]; 4879371c9d4SSatish Balay } else if (cone[nc - 3] == points[0]) { 4889371c9d4SSatish Balay tmp = cone[nc - 3]; 4899371c9d4SSatish Balay cone[nc - 3] = cone[nc - 1]; 4909371c9d4SSatish Balay cone[nc - 1] = tmp; 4919371c9d4SSatish Balay cone[nc++] = points[1]; 4929371c9d4SSatish Balay if (cone[0] != points[2]) cone[nc++] = points[2]; 4939371c9d4SSatish Balay } else if (cone[nc - 3] == points[2]) { 4949371c9d4SSatish Balay tmp = cone[nc - 3]; 4959371c9d4SSatish Balay cone[nc - 3] = cone[nc - 1]; 4969371c9d4SSatish Balay cone[nc - 1] = tmp; 4979371c9d4SSatish Balay cone[nc++] = points[1]; 4989371c9d4SSatish Balay if (cone[0] != points[0]) cone[nc++] = points[0]; 4999371c9d4SSatish Balay } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %d does not match its predecessor", eid); 5007bee2925SMatthew Knepley } 5017bee2925SMatthew Knepley } 50263a3b9bcSJacob Faibussowitsch PetscCheck(nc == 2 * Ner, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %" PetscInt_FMT " != %" PetscInt_FMT, nc, 2 * Ner); 503ad540459SPierre Jolivet if (Ner == 4) cone[nc++] = numVertices - newVertices + numEdges + numQuads++; 50463a3b9bcSJacob Faibussowitsch PetscCheck(nc <= maxCorners, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %" PetscInt_FMT " > %" PetscInt_FMT " max", nc, maxCorners); 5057bee2925SMatthew Knepley /* Triangulate the loop */ 5067bee2925SMatthew Knepley switch (Ner) { 5077bee2925SMatthew Knepley case 2: /* Bi-Segment -> 2 triangles */ 5087bee2925SMatthew Knepley Nt = 2; 5097bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[0]; 5107bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[1]; 5117bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[2]; 5127bee2925SMatthew Knepley ++cOff; 5137bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[0]; 5147bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[2]; 5157bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[3]; 5167bee2925SMatthew Knepley ++cOff; 5177bee2925SMatthew Knepley break; 5187bee2925SMatthew Knepley case 3: /* Triangle -> 4 triangles */ 5197bee2925SMatthew Knepley Nt = 4; 5207bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[0]; 5217bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[1]; 5227bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[5]; 5237bee2925SMatthew Knepley ++cOff; 5247bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[1]; 5257bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[2]; 5267bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[3]; 5277bee2925SMatthew Knepley ++cOff; 5287bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[5]; 5297bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[3]; 5307bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[4]; 5317bee2925SMatthew Knepley ++cOff; 5327bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[1]; 5337bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[3]; 5347bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[5]; 5357bee2925SMatthew Knepley ++cOff; 5367bee2925SMatthew Knepley break; 5377bee2925SMatthew Knepley case 4: /* Quad -> 8 triangles */ 5387bee2925SMatthew Knepley Nt = 8; 5397bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[0]; 5407bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[1]; 5417bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[7]; 5427bee2925SMatthew Knepley ++cOff; 5437bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[1]; 5447bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[2]; 5457bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[3]; 5467bee2925SMatthew Knepley ++cOff; 5477bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[3]; 5487bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[4]; 5497bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[5]; 5507bee2925SMatthew Knepley ++cOff; 5517bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[5]; 5527bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[6]; 5537bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[7]; 5547bee2925SMatthew Knepley ++cOff; 5557bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[8]; 5567bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[1]; 5577bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[3]; 5587bee2925SMatthew Knepley ++cOff; 5597bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[8]; 5607bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[3]; 5617bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[5]; 5627bee2925SMatthew Knepley ++cOff; 5637bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[8]; 5647bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[5]; 5657bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[7]; 5667bee2925SMatthew Knepley ++cOff; 5677bee2925SMatthew Knepley cells[cOff * numCorners + 0] = cone[8]; 5687bee2925SMatthew Knepley cells[cOff * numCorners + 1] = cone[7]; 5697bee2925SMatthew Knepley cells[cOff * numCorners + 2] = cone[1]; 5707bee2925SMatthew Knepley ++cOff; 5717bee2925SMatthew Knepley break; 572d71ae5a4SJacob Faibussowitsch default: 573d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d edges, which we do not support", lid, Ner); 5747bee2925SMatthew Knepley } 575266cfabeSMatthew G. Knepley if (debug) { 5767bee2925SMatthew Knepley for (t = 0; t < Nt; ++t) { 57763a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " LOOP Corner NODEs Triangle %" PetscInt_FMT " (", t)); 5787bee2925SMatthew Knepley for (c = 0; c < numCorners; ++c) { 5799566063dSJacob Faibussowitsch if (c > 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, ", ")); 58063a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "%" PetscInt_FMT, cells[(cOff - Nt + t) * numCorners + c])); 5817bee2925SMatthew Knepley } 5829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, ")\n")); 5837bee2925SMatthew Knepley } 5847bee2925SMatthew Knepley } 585266cfabeSMatthew G. Knepley } 5867bee2925SMatthew Knepley EG_free(lobjs); 5877bee2925SMatthew Knepley } 5887bee2925SMatthew Knepley } 58963a3b9bcSJacob Faibussowitsch PetscCheck(cOff == numCells, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count of total cells %" PetscInt_FMT " != %" PetscInt_FMT " previous count", cOff, numCells); 5909566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numVertices, numCorners, PETSC_TRUE, cells, cdim, coords, &dm)); 5919566063dSJacob Faibussowitsch PetscCall(PetscFree3(coords, cells, cone)); 59263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(dm, " Total Number of Unique Cells = %" PetscInt_FMT " (%" PetscInt_FMT ")\n", numCells, newCells)); 59363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(dm, " Total Number of Unique Vertices = %" PetscInt_FMT " (%" PetscInt_FMT ")\n", numVertices, newVertices)); 5947bee2925SMatthew Knepley /* Embed EGADS model in DM */ 5957bee2925SMatthew Knepley { 5967bee2925SMatthew Knepley PetscContainer modelObj, contextObj; 5977bee2925SMatthew Knepley 5989566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &modelObj)); 5999566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(modelObj, model)); 6009566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Model", (PetscObject)modelObj)); 6019566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&modelObj)); 6027bee2925SMatthew Knepley 6039566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &contextObj)); 6049566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(contextObj, context)); 6059566063dSJacob Faibussowitsch PetscCall(PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private)); 6069566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Context", (PetscObject)contextObj)); 6079566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&contextObj)); 6087bee2925SMatthew Knepley } 6097bee2925SMatthew Knepley /* Label points */ 6109566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Body ID")); 6119566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 6129566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Face ID")); 6139566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 6149566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Edge ID")); 6159566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel)); 6169566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Vertex ID")); 6179566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel)); 6187bee2925SMatthew Knepley cOff = 0; 6197bee2925SMatthew Knepley for (b = 0; b < nbodies; ++b) { 6207bee2925SMatthew Knepley ego body = bodies[b]; 6217bee2925SMatthew Knepley int id, Nl, l; 6227bee2925SMatthew Knepley 6239566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs)); 6247bee2925SMatthew Knepley for (l = 0; l < Nl; ++l) { 6257bee2925SMatthew Knepley ego loop = lobjs[l]; 6267bee2925SMatthew Knepley ego *fobjs; 6277bee2925SMatthew Knepley int lid, Nf, fid, Ner = 0, Ne, e, Nt = 0, t; 6287bee2925SMatthew Knepley 6295f80ce2aSJacob Faibussowitsch lid = EG_indexBodyTopo(body, loop); 6309566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs)); 63108401ef6SPierre Jolivet PetscCheck(Nf <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf); 6325f80ce2aSJacob Faibussowitsch fid = EG_indexBodyTopo(body, fobjs[0]); 6337bee2925SMatthew Knepley EG_free(fobjs); 6349566063dSJacob Faibussowitsch PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses)); 6357bee2925SMatthew Knepley for (e = 0; e < Ne; ++e) { 6367bee2925SMatthew Knepley ego edge = objs[e]; 6377bee2925SMatthew Knepley int eid, Nv, v; 6387bee2925SMatthew Knepley PetscInt points[3], support[2], numEdges, edgeNum; 6397bee2925SMatthew Knepley const PetscInt *edges; 6407bee2925SMatthew Knepley 6417bee2925SMatthew Knepley eid = EG_indexBodyTopo(body, edge); 6429566063dSJacob Faibussowitsch PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 6437bee2925SMatthew Knepley if (mtype == DEGENERATE) continue; 6447bee2925SMatthew Knepley else ++Ner; 6457bee2925SMatthew Knepley for (v = 0; v < Nv; ++v) { 6467bee2925SMatthew Knepley ego vertex = nobjs[v]; 6477bee2925SMatthew Knepley 6487bee2925SMatthew Knepley id = EG_indexBodyTopo(body, vertex); 6499566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(edgeLabel, numCells + id - 1, eid)); 6507bee2925SMatthew Knepley points[v * 2] = numCells + id - 1; 6517bee2925SMatthew Knepley } 6529566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(edgeMap, eid - 1, &edgeNum)); 6537bee2925SMatthew Knepley points[1] = numCells + numVertices - newVertices + edgeNum; 6547bee2925SMatthew Knepley 6559566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(edgeLabel, points[1], eid)); 6567bee2925SMatthew Knepley support[0] = points[0]; 6577bee2925SMatthew Knepley support[1] = points[1]; 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 support[0] = points[1]; 6637bee2925SMatthew Knepley support[1] = points[2]; 6649566063dSJacob Faibussowitsch PetscCall(DMPlexGetJoin(dm, 2, support, &numEdges, &edges)); 66563a3b9bcSJacob 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); 6669566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(edgeLabel, edges[0], eid)); 6679566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges)); 6687bee2925SMatthew Knepley } 6697bee2925SMatthew Knepley switch (Ner) { 670d71ae5a4SJacob Faibussowitsch case 2: 671d71ae5a4SJacob Faibussowitsch Nt = 2; 672d71ae5a4SJacob Faibussowitsch break; 673d71ae5a4SJacob Faibussowitsch case 3: 674d71ae5a4SJacob Faibussowitsch Nt = 4; 675d71ae5a4SJacob Faibussowitsch break; 676d71ae5a4SJacob Faibussowitsch case 4: 677d71ae5a4SJacob Faibussowitsch Nt = 8; 678d71ae5a4SJacob Faibussowitsch break; 679d71ae5a4SJacob Faibussowitsch default: 680d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Loop with %d edges is unsupported", Ner); 6817bee2925SMatthew Knepley } 6827bee2925SMatthew Knepley for (t = 0; t < Nt; ++t) { 6839566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, cOff + t, b)); 6849566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(faceLabel, cOff + t, fid)); 6857bee2925SMatthew Knepley } 6867bee2925SMatthew Knepley cOff += Nt; 6877bee2925SMatthew Knepley } 6887bee2925SMatthew Knepley EG_free(lobjs); 6897bee2925SMatthew Knepley } 6909566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&edgeMap)); 6919566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 6927bee2925SMatthew Knepley for (c = cStart; c < cEnd; ++c) { 6937bee2925SMatthew Knepley PetscInt *closure = NULL; 6947bee2925SMatthew Knepley PetscInt clSize, cl, bval, fval; 6957bee2925SMatthew Knepley 6969566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure)); 6979566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, c, &bval)); 6989566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(faceLabel, c, &fval)); 6997bee2925SMatthew Knepley for (cl = 0; cl < clSize * 2; cl += 2) { 7009566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, closure[cl], bval)); 7019566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(faceLabel, closure[cl], fval)); 7027bee2925SMatthew Knepley } 7039566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure)); 7047bee2925SMatthew Knepley } 7057bee2925SMatthew Knepley *newdm = dm; 7067bee2925SMatthew Knepley PetscFunctionReturn(0); 7077bee2925SMatthew Knepley } 708c1cad2e7SMatthew G. Knepley 709d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateEGADS(MPI_Comm comm, ego context, ego model, DM *newdm) 710d71ae5a4SJacob Faibussowitsch { 711c1cad2e7SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel; 712c1cad2e7SMatthew G. Knepley // EGADS/EGADSLite variables 713c1cad2e7SMatthew G. Knepley ego geom, *bodies, *mobjs, *fobjs, *lobjs, *eobjs, *nobjs; 714c1cad2e7SMatthew G. Knepley ego topRef, prev, next; 715c1cad2e7SMatthew G. Knepley int oclass, mtype, nbodies, *senses, *lSenses, *eSenses; 716c1cad2e7SMatthew G. Knepley int b; 717c1cad2e7SMatthew G. Knepley // PETSc variables 718c1cad2e7SMatthew G. Knepley DM dm; 719c1cad2e7SMatthew G. Knepley PetscHMapI edgeMap = NULL, bodyIndexMap = NULL, bodyVertexMap = NULL, bodyEdgeMap = NULL, bodyFaceMap = NULL, bodyEdgeGlobalMap = NULL; 720c1cad2e7SMatthew G. Knepley PetscInt dim = -1, cdim = -1, numCorners = 0, numVertices = 0, numEdges = 0, numFaces = 0, numCells = 0, edgeCntr = 0; 721c1cad2e7SMatthew G. Knepley PetscInt cellCntr = 0, numPoints = 0; 722c1cad2e7SMatthew G. Knepley PetscInt *cells = NULL; 723c1cad2e7SMatthew G. Knepley const PetscInt *cone = NULL; 724c1cad2e7SMatthew G. Knepley PetscReal *coords = NULL; 725c1cad2e7SMatthew G. Knepley PetscMPIInt rank; 726c1cad2e7SMatthew G. Knepley 727c1cad2e7SMatthew G. Knepley PetscFunctionBeginUser; 7289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 729c5853193SPierre Jolivet if (rank == 0) { 730c1cad2e7SMatthew G. Knepley // --------------------------------------------------------------------------------------------------- 731c1cad2e7SMatthew G. Knepley // Generate Petsc Plex 732c1cad2e7SMatthew G. Knepley // Get all Nodes in model, record coordinates in a correctly formatted array 733c1cad2e7SMatthew G. Knepley // Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array 734c1cad2e7SMatthew G. Knepley // We need to uniformly refine the initial geometry to guarantee a valid mesh 735c1cad2e7SMatthew G. Knepley 736*d5b43468SJose E. Roman // Calculate cell and vertex sizes 7379566063dSJacob Faibussowitsch PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses)); 738c1cad2e7SMatthew G. Knepley 7399566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&edgeMap)); 7409566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&bodyIndexMap)); 7419566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&bodyVertexMap)); 7429566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&bodyEdgeMap)); 7439566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&bodyEdgeGlobalMap)); 7449566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&bodyFaceMap)); 745c1cad2e7SMatthew G. Knepley 746c1cad2e7SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 747c1cad2e7SMatthew G. Knepley ego body = bodies[b]; 748c1cad2e7SMatthew G. Knepley int Nf, Ne, Nv; 749c1cad2e7SMatthew G. Knepley PetscHashIter BIiter, BViter, BEiter, BEGiter, BFiter, EMiter; 750c1cad2e7SMatthew G. Knepley PetscBool BIfound, BVfound, BEfound, BEGfound, BFfound, EMfound; 751c1cad2e7SMatthew G. Knepley 7529566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyIndexMap, b, &BIiter, &BIfound)); 7539566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound)); 7549566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound)); 7559566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound)); 7569566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound)); 757c1cad2e7SMatthew G. Knepley 7589566063dSJacob Faibussowitsch if (!BIfound) PetscCall(PetscHMapISet(bodyIndexMap, b, numFaces + numEdges + numVertices)); 7599566063dSJacob Faibussowitsch if (!BVfound) PetscCall(PetscHMapISet(bodyVertexMap, b, numVertices)); 7609566063dSJacob Faibussowitsch if (!BEfound) PetscCall(PetscHMapISet(bodyEdgeMap, b, numEdges)); 7619566063dSJacob Faibussowitsch if (!BEGfound) PetscCall(PetscHMapISet(bodyEdgeGlobalMap, b, edgeCntr)); 7629566063dSJacob Faibussowitsch if (!BFfound) PetscCall(PetscHMapISet(bodyFaceMap, b, numFaces)); 763c1cad2e7SMatthew G. Knepley 7649566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 7659566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs)); 7669566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs)); 767c1cad2e7SMatthew G. Knepley EG_free(fobjs); 768c1cad2e7SMatthew G. Knepley EG_free(eobjs); 769c1cad2e7SMatthew G. Knepley EG_free(nobjs); 770c1cad2e7SMatthew G. Knepley 771c1cad2e7SMatthew G. Knepley // Remove DEGENERATE EDGES from Edge count 7729566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs)); 773c1cad2e7SMatthew G. Knepley int Netemp = 0; 774c1cad2e7SMatthew G. Knepley for (int e = 0; e < Ne; ++e) { 775c1cad2e7SMatthew G. Knepley ego edge = eobjs[e]; 776c1cad2e7SMatthew G. Knepley int eid; 777c1cad2e7SMatthew G. Knepley 7789566063dSJacob Faibussowitsch PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next)); 7795f80ce2aSJacob Faibussowitsch eid = EG_indexBodyTopo(body, edge); 780c1cad2e7SMatthew G. Knepley 7819566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(edgeMap, edgeCntr + eid - 1, &EMiter, &EMfound)); 782c1cad2e7SMatthew G. Knepley if (mtype == DEGENERATE) { 7839566063dSJacob Faibussowitsch if (!EMfound) PetscCall(PetscHMapISet(edgeMap, edgeCntr + eid - 1, -1)); 7849371c9d4SSatish Balay } else { 785c1cad2e7SMatthew G. Knepley ++Netemp; 7869566063dSJacob Faibussowitsch if (!EMfound) PetscCall(PetscHMapISet(edgeMap, edgeCntr + eid - 1, Netemp)); 787c1cad2e7SMatthew G. Knepley } 788c1cad2e7SMatthew G. Knepley } 789c1cad2e7SMatthew G. Knepley EG_free(eobjs); 790c1cad2e7SMatthew G. Knepley 791c1cad2e7SMatthew G. Knepley // Determine Number of Cells 7929566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 793c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 794c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 795c1cad2e7SMatthew G. Knepley int edgeTemp = 0; 796c1cad2e7SMatthew G. Knepley 7979566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, face, EDGE, &Ne, &eobjs)); 798c1cad2e7SMatthew G. Knepley for (int e = 0; e < Ne; ++e) { 799c1cad2e7SMatthew G. Knepley ego edge = eobjs[e]; 800c1cad2e7SMatthew G. Knepley 8019566063dSJacob Faibussowitsch PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next)); 802ad540459SPierre Jolivet if (mtype != DEGENERATE) ++edgeTemp; 803c1cad2e7SMatthew G. Knepley } 804c1cad2e7SMatthew G. Knepley numCells += (2 * edgeTemp); 805c1cad2e7SMatthew G. Knepley EG_free(eobjs); 806c1cad2e7SMatthew G. Knepley } 807c1cad2e7SMatthew G. Knepley EG_free(fobjs); 808c1cad2e7SMatthew G. Knepley 809c1cad2e7SMatthew G. Knepley numFaces += Nf; 810c1cad2e7SMatthew G. Knepley numEdges += Netemp; 811c1cad2e7SMatthew G. Knepley numVertices += Nv; 812c1cad2e7SMatthew G. Knepley edgeCntr += Ne; 813c1cad2e7SMatthew G. Knepley } 814c1cad2e7SMatthew G. Knepley 815c1cad2e7SMatthew G. Knepley // Set up basic DMPlex parameters 816c1cad2e7SMatthew G. Knepley dim = 2; // Assumes 3D Models :: Need to handle 2D modles in the future 817c1cad2e7SMatthew G. Knepley cdim = 3; // Assumes 3D Models :: Need to update to handle 2D modles in future 818c1cad2e7SMatthew G. Knepley numCorners = 3; // Split Faces into triangles 819c1cad2e7SMatthew G. Knepley numPoints = numVertices + numEdges + numFaces; // total number of coordinate points 820c1cad2e7SMatthew G. Knepley 8219566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numPoints * cdim, &coords, numCells * numCorners, &cells)); 822c1cad2e7SMatthew G. Knepley 823c1cad2e7SMatthew G. Knepley // Get Vertex Coordinates and Set up Cells 824c1cad2e7SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 825c1cad2e7SMatthew G. Knepley ego body = bodies[b]; 826c1cad2e7SMatthew G. Knepley int Nf, Ne, Nv; 827c1cad2e7SMatthew G. Knepley PetscInt bodyVertexIndexStart, bodyEdgeIndexStart, bodyEdgeGlobalIndexStart, bodyFaceIndexStart; 828c1cad2e7SMatthew G. Knepley PetscHashIter BViter, BEiter, BEGiter, BFiter, EMiter; 829c1cad2e7SMatthew G. Knepley PetscBool BVfound, BEfound, BEGfound, BFfound, EMfound; 830c1cad2e7SMatthew G. Knepley 831c1cad2e7SMatthew G. Knepley // Vertices on Current Body 8329566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs)); 833c1cad2e7SMatthew G. Knepley 8349566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound)); 83528b400f6SJacob Faibussowitsch PetscCheck(BVfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyVertexMap", b); 8369566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(bodyVertexMap, b, &bodyVertexIndexStart)); 837c1cad2e7SMatthew G. Knepley 838c1cad2e7SMatthew G. Knepley for (int v = 0; v < Nv; ++v) { 839c1cad2e7SMatthew G. Knepley ego vertex = nobjs[v]; 840c1cad2e7SMatthew G. Knepley double limits[4]; 841c1cad2e7SMatthew G. Knepley int id, dummy; 842c1cad2e7SMatthew G. Knepley 8439566063dSJacob Faibussowitsch PetscCall(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses)); 8445f80ce2aSJacob Faibussowitsch id = EG_indexBodyTopo(body, vertex); 845c1cad2e7SMatthew G. Knepley 846c1cad2e7SMatthew G. Knepley coords[(bodyVertexIndexStart + id - 1) * cdim + 0] = limits[0]; 847c1cad2e7SMatthew G. Knepley coords[(bodyVertexIndexStart + id - 1) * cdim + 1] = limits[1]; 848c1cad2e7SMatthew G. Knepley coords[(bodyVertexIndexStart + id - 1) * cdim + 2] = limits[2]; 849c1cad2e7SMatthew G. Knepley } 850c1cad2e7SMatthew G. Knepley EG_free(nobjs); 851c1cad2e7SMatthew G. Knepley 852c1cad2e7SMatthew G. Knepley // Edge Midpoint Vertices on Current Body 8539566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs)); 854c1cad2e7SMatthew G. Knepley 8559566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound)); 85628b400f6SJacob Faibussowitsch PetscCheck(BEfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeMap", b); 8579566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(bodyEdgeMap, b, &bodyEdgeIndexStart)); 858c1cad2e7SMatthew G. Knepley 8599566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound)); 86028b400f6SJacob Faibussowitsch PetscCheck(BEGfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeGlobalMap", b); 8619566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(bodyEdgeGlobalMap, b, &bodyEdgeGlobalIndexStart)); 862c1cad2e7SMatthew G. Knepley 863c1cad2e7SMatthew G. Knepley for (int e = 0; e < Ne; ++e) { 864c1cad2e7SMatthew G. Knepley ego edge = eobjs[e]; 865c1cad2e7SMatthew G. Knepley double range[2], avgt[1], cntrPnt[9]; 866c1cad2e7SMatthew G. Knepley int eid, eOffset; 867c1cad2e7SMatthew G. Knepley int periodic; 868c1cad2e7SMatthew G. Knepley 8699566063dSJacob Faibussowitsch PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next)); 870ad540459SPierre Jolivet if (mtype == DEGENERATE) continue; 871c1cad2e7SMatthew G. Knepley 8725f80ce2aSJacob Faibussowitsch eid = EG_indexBodyTopo(body, edge); 873c1cad2e7SMatthew G. Knepley 874c1cad2e7SMatthew G. Knepley // get relative offset from globalEdgeID Vector 8759566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound)); 87628b400f6SJacob Faibussowitsch PetscCheck(EMfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %d not found in edgeMap", bodyEdgeGlobalIndexStart + eid - 1); 8779566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset)); 878c1cad2e7SMatthew G. Knepley 8799566063dSJacob Faibussowitsch PetscCall(EG_getRange(edge, range, &periodic)); 880c1cad2e7SMatthew G. Knepley avgt[0] = (range[0] + range[1]) / 2.; 881c1cad2e7SMatthew G. Knepley 8829566063dSJacob Faibussowitsch PetscCall(EG_evaluate(edge, avgt, cntrPnt)); 883c1cad2e7SMatthew G. Knepley coords[(numVertices + bodyEdgeIndexStart + eOffset - 1) * cdim + 0] = cntrPnt[0]; 884c1cad2e7SMatthew G. Knepley coords[(numVertices + bodyEdgeIndexStart + eOffset - 1) * cdim + 1] = cntrPnt[1]; 885c1cad2e7SMatthew G. Knepley coords[(numVertices + bodyEdgeIndexStart + eOffset - 1) * cdim + 2] = cntrPnt[2]; 886c1cad2e7SMatthew G. Knepley } 887c1cad2e7SMatthew G. Knepley EG_free(eobjs); 888c1cad2e7SMatthew G. Knepley 889c1cad2e7SMatthew G. Knepley // Face Midpoint Vertices on Current Body 8909566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 891c1cad2e7SMatthew G. Knepley 8929566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound)); 89328b400f6SJacob Faibussowitsch PetscCheck(BFfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyFaceMap", b); 8949566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(bodyFaceMap, b, &bodyFaceIndexStart)); 895c1cad2e7SMatthew G. Knepley 896c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 897c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 898c1cad2e7SMatthew G. Knepley double range[4], avgUV[2], cntrPnt[18]; 899c1cad2e7SMatthew G. Knepley int peri, id; 900c1cad2e7SMatthew G. Knepley 901c1cad2e7SMatthew G. Knepley id = EG_indexBodyTopo(body, face); 9029566063dSJacob Faibussowitsch PetscCall(EG_getRange(face, range, &peri)); 903c1cad2e7SMatthew G. Knepley 904c1cad2e7SMatthew G. Knepley avgUV[0] = (range[0] + range[1]) / 2.; 905c1cad2e7SMatthew G. Knepley avgUV[1] = (range[2] + range[3]) / 2.; 9069566063dSJacob Faibussowitsch PetscCall(EG_evaluate(face, avgUV, cntrPnt)); 907c1cad2e7SMatthew G. Knepley 908c1cad2e7SMatthew G. Knepley coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1) * cdim + 0] = cntrPnt[0]; 909c1cad2e7SMatthew G. Knepley coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1) * cdim + 1] = cntrPnt[1]; 910c1cad2e7SMatthew G. Knepley coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1) * cdim + 2] = cntrPnt[2]; 911c1cad2e7SMatthew G. Knepley } 912c1cad2e7SMatthew G. Knepley EG_free(fobjs); 913c1cad2e7SMatthew G. Knepley 914c1cad2e7SMatthew G. Knepley // Define Cells :: Note - This could be incorporated in the Face Midpoint Vertices Loop but was kept separate for clarity 9159566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 916c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 917c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 918c1cad2e7SMatthew G. Knepley int fID, midFaceID, midPntID, startID, endID, Nl; 919c1cad2e7SMatthew G. Knepley 9205f80ce2aSJacob Faibussowitsch fID = EG_indexBodyTopo(body, face); 921c1cad2e7SMatthew G. Knepley midFaceID = numVertices + numEdges + bodyFaceIndexStart + fID - 1; 922c1cad2e7SMatthew G. Knepley // Must Traverse Loop to ensure we have all necessary information like the sense (+/- 1) of the edges. 923c1cad2e7SMatthew G. Knepley // TODO :: Only handles single loop faces (No holes). The choices for handling multiloop faces are: 924c1cad2e7SMatthew G. Knepley // 1) Use the DMPlexCreateEGADSFromFile() with the -dm_plex_egads_with_tess = 1 option. 925c1cad2e7SMatthew G. Knepley // This will use a default EGADS tessellation as an initial surface mesh. 926*d5b43468SJose E. Roman // 2) Create the initial surface mesh via a 2D mesher :: Currently not available (?future?) 927c1cad2e7SMatthew G. Knepley // May I suggest the XXXX as a starting point? 928c1cad2e7SMatthew G. Knepley 9299566063dSJacob Faibussowitsch PetscCall(EG_getTopology(face, &geom, &oclass, &mtype, NULL, &Nl, &lobjs, &lSenses)); 930c1cad2e7SMatthew G. Knepley 93108401ef6SPierre 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); 932c1cad2e7SMatthew G. Knepley for (int l = 0; l < Nl; ++l) { 933c1cad2e7SMatthew G. Knepley ego loop = lobjs[l]; 934c1cad2e7SMatthew G. Knepley 9359566063dSJacob Faibussowitsch PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses)); 936c1cad2e7SMatthew G. Knepley for (int e = 0; e < Ne; ++e) { 937c1cad2e7SMatthew G. Knepley ego edge = eobjs[e]; 938c1cad2e7SMatthew G. Knepley int eid, eOffset; 939c1cad2e7SMatthew G. Knepley 9409566063dSJacob Faibussowitsch PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next)); 941c1cad2e7SMatthew G. Knepley eid = EG_indexBodyTopo(body, edge); 942ad540459SPierre Jolivet if (mtype == DEGENERATE) continue; 943c1cad2e7SMatthew G. Knepley 944c1cad2e7SMatthew G. Knepley // get relative offset from globalEdgeID Vector 9459566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound)); 94628b400f6SJacob 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); 9479566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset)); 948c1cad2e7SMatthew G. Knepley 949c1cad2e7SMatthew G. Knepley midPntID = numVertices + bodyEdgeIndexStart + eOffset - 1; 950c1cad2e7SMatthew G. Knepley 9519566063dSJacob Faibussowitsch PetscCall(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 952c1cad2e7SMatthew G. Knepley 9539371c9d4SSatish Balay if (eSenses[e] > 0) { 9549371c9d4SSatish Balay startID = EG_indexBodyTopo(body, nobjs[0]); 9559371c9d4SSatish Balay endID = EG_indexBodyTopo(body, nobjs[1]); 9569371c9d4SSatish Balay } else { 9579371c9d4SSatish Balay startID = EG_indexBodyTopo(body, nobjs[1]); 9589371c9d4SSatish Balay endID = EG_indexBodyTopo(body, nobjs[0]); 9599371c9d4SSatish Balay } 960c1cad2e7SMatthew G. Knepley 961c1cad2e7SMatthew G. Knepley // Define 2 Cells per Edge with correct orientation 962c1cad2e7SMatthew G. Knepley cells[cellCntr * numCorners + 0] = midFaceID; 963c1cad2e7SMatthew G. Knepley cells[cellCntr * numCorners + 1] = bodyVertexIndexStart + startID - 1; 964c1cad2e7SMatthew G. Knepley cells[cellCntr * numCorners + 2] = midPntID; 965c1cad2e7SMatthew G. Knepley 966c1cad2e7SMatthew G. Knepley cells[cellCntr * numCorners + 3] = midFaceID; 967c1cad2e7SMatthew G. Knepley cells[cellCntr * numCorners + 4] = midPntID; 968c1cad2e7SMatthew G. Knepley cells[cellCntr * numCorners + 5] = bodyVertexIndexStart + endID - 1; 969c1cad2e7SMatthew G. Knepley 970c1cad2e7SMatthew G. Knepley cellCntr = cellCntr + 2; 971c1cad2e7SMatthew G. Knepley } 972c1cad2e7SMatthew G. Knepley } 973c1cad2e7SMatthew G. Knepley } 974c1cad2e7SMatthew G. Knepley EG_free(fobjs); 975c1cad2e7SMatthew G. Knepley } 976c1cad2e7SMatthew G. Knepley } 977c1cad2e7SMatthew G. Knepley 978c1cad2e7SMatthew G. Knepley // Generate DMPlex 9799566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numPoints, numCorners, PETSC_TRUE, cells, cdim, coords, &dm)); 9809566063dSJacob Faibussowitsch PetscCall(PetscFree2(coords, cells)); 98163a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(dm, " Total Number of Unique Cells = %" PetscInt_FMT " \n", numCells)); 98263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(dm, " Total Number of Unique Vertices = %" PetscInt_FMT " \n", numVertices)); 983c1cad2e7SMatthew G. Knepley 984c1cad2e7SMatthew G. Knepley // Embed EGADS model in DM 985c1cad2e7SMatthew G. Knepley { 986c1cad2e7SMatthew G. Knepley PetscContainer modelObj, contextObj; 987c1cad2e7SMatthew G. Knepley 9889566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &modelObj)); 9899566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(modelObj, model)); 9909566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Model", (PetscObject)modelObj)); 9919566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&modelObj)); 992c1cad2e7SMatthew G. Knepley 9939566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &contextObj)); 9949566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(contextObj, context)); 9959566063dSJacob Faibussowitsch PetscCall(PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private)); 9969566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Context", (PetscObject)contextObj)); 9979566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&contextObj)); 998c1cad2e7SMatthew G. Knepley } 999c1cad2e7SMatthew G. Knepley // Label points 1000c1cad2e7SMatthew G. Knepley PetscInt nStart, nEnd; 1001c1cad2e7SMatthew G. Knepley 10029566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Body ID")); 10039566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 10049566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Face ID")); 10059566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 10069566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Edge ID")); 10079566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel)); 10089566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Vertex ID")); 10099566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel)); 1010c1cad2e7SMatthew G. Knepley 10119566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 2, &nStart, &nEnd)); 1012c1cad2e7SMatthew G. Knepley 1013c1cad2e7SMatthew G. Knepley cellCntr = 0; 1014c1cad2e7SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 1015c1cad2e7SMatthew G. Knepley ego body = bodies[b]; 1016c1cad2e7SMatthew G. Knepley int Nv, Ne, Nf; 1017c1cad2e7SMatthew G. Knepley PetscInt bodyVertexIndexStart, bodyEdgeIndexStart, bodyEdgeGlobalIndexStart, bodyFaceIndexStart; 1018c1cad2e7SMatthew G. Knepley PetscHashIter BViter, BEiter, BEGiter, BFiter, EMiter; 1019c1cad2e7SMatthew G. Knepley PetscBool BVfound, BEfound, BEGfound, BFfound, EMfound; 1020c1cad2e7SMatthew G. Knepley 10219566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound)); 102228b400f6SJacob Faibussowitsch PetscCheck(BVfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyVertexMap", b); 10239566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(bodyVertexMap, b, &bodyVertexIndexStart)); 1024c1cad2e7SMatthew G. Knepley 10259566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound)); 102628b400f6SJacob Faibussowitsch PetscCheck(BEfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeMap", b); 10279566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(bodyEdgeMap, b, &bodyEdgeIndexStart)); 1028c1cad2e7SMatthew G. Knepley 10299566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound)); 103028b400f6SJacob Faibussowitsch PetscCheck(BFfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyFaceMap", b); 10319566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(bodyFaceMap, b, &bodyFaceIndexStart)); 1032c1cad2e7SMatthew G. Knepley 10339566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound)); 103428b400f6SJacob Faibussowitsch PetscCheck(BEGfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeGlobalMap", b); 10359566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(bodyEdgeGlobalMap, b, &bodyEdgeGlobalIndexStart)); 1036c1cad2e7SMatthew G. Knepley 10379566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 1038c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 1039c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 1040c1cad2e7SMatthew G. Knepley int fID, Nl; 1041c1cad2e7SMatthew G. Knepley 10425f80ce2aSJacob Faibussowitsch fID = EG_indexBodyTopo(body, face); 1043c1cad2e7SMatthew G. Knepley 10449566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, face, LOOP, &Nl, &lobjs)); 1045c1cad2e7SMatthew G. Knepley for (int l = 0; l < Nl; ++l) { 1046c1cad2e7SMatthew G. Knepley ego loop = lobjs[l]; 1047c1cad2e7SMatthew G. Knepley int lid; 1048c1cad2e7SMatthew G. Knepley 10495f80ce2aSJacob Faibussowitsch lid = EG_indexBodyTopo(body, loop); 105008401ef6SPierre Jolivet PetscCheck(Nl <= 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf); 1051c1cad2e7SMatthew G. Knepley 10529566063dSJacob Faibussowitsch PetscCall(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses)); 1053c1cad2e7SMatthew G. Knepley for (int e = 0; e < Ne; ++e) { 1054c1cad2e7SMatthew G. Knepley ego edge = eobjs[e]; 1055c1cad2e7SMatthew G. Knepley int eid, eOffset; 1056c1cad2e7SMatthew G. Knepley 1057c1cad2e7SMatthew G. Knepley // Skip DEGENERATE Edges 10589566063dSJacob Faibussowitsch PetscCall(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next)); 1059ad540459SPierre Jolivet if (mtype == DEGENERATE) continue; 10605f80ce2aSJacob Faibussowitsch eid = EG_indexBodyTopo(body, edge); 1061c1cad2e7SMatthew G. Knepley 1062c1cad2e7SMatthew G. Knepley // get relative offset from globalEdgeID Vector 10639566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound)); 106428b400f6SJacob 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); 10659566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset)); 1066c1cad2e7SMatthew G. Knepley 10679566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, edge, NODE, &Nv, &nobjs)); 1068c1cad2e7SMatthew G. Knepley for (int v = 0; v < Nv; ++v) { 1069c1cad2e7SMatthew G. Knepley ego vertex = nobjs[v]; 1070c1cad2e7SMatthew G. Knepley int vID; 1071c1cad2e7SMatthew G. Knepley 10725f80ce2aSJacob Faibussowitsch vID = EG_indexBodyTopo(body, vertex); 10739566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, nStart + bodyVertexIndexStart + vID - 1, b)); 10749566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(vertexLabel, nStart + bodyVertexIndexStart + vID - 1, vID)); 1075c1cad2e7SMatthew G. Knepley } 1076c1cad2e7SMatthew G. Knepley EG_free(nobjs); 1077c1cad2e7SMatthew G. Knepley 10789566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, nStart + numVertices + bodyEdgeIndexStart + eOffset - 1, b)); 10799566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(edgeLabel, nStart + numVertices + bodyEdgeIndexStart + eOffset - 1, eid)); 1080c1cad2e7SMatthew G. Knepley 1081c1cad2e7SMatthew G. Knepley // Define Cell faces 1082c1cad2e7SMatthew G. Knepley for (int jj = 0; jj < 2; ++jj) { 10839566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, cellCntr, b)); 10849566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(faceLabel, cellCntr, fID)); 10859566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, cellCntr, &cone)); 1086c1cad2e7SMatthew G. Knepley 10879566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, cone[0], b)); 10889566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(faceLabel, cone[0], fID)); 1089c1cad2e7SMatthew G. Knepley 10909566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, cone[1], b)); 10919566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(edgeLabel, cone[1], eid)); 1092c1cad2e7SMatthew G. Knepley 10939566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, cone[2], b)); 10949566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(faceLabel, cone[2], fID)); 1095c1cad2e7SMatthew G. Knepley 1096c1cad2e7SMatthew G. Knepley cellCntr = cellCntr + 1; 1097c1cad2e7SMatthew G. Knepley } 1098c1cad2e7SMatthew G. Knepley } 1099c1cad2e7SMatthew G. Knepley } 1100c1cad2e7SMatthew G. Knepley EG_free(lobjs); 1101c1cad2e7SMatthew G. Knepley 11029566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, nStart + numVertices + numEdges + bodyFaceIndexStart + fID - 1, b)); 11039566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(faceLabel, nStart + numVertices + numEdges + bodyFaceIndexStart + fID - 1, fID)); 1104c1cad2e7SMatthew G. Knepley } 1105c1cad2e7SMatthew G. Knepley EG_free(fobjs); 1106c1cad2e7SMatthew G. Knepley } 1107c1cad2e7SMatthew G. Knepley 11089566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&edgeMap)); 11099566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&bodyIndexMap)); 11109566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&bodyVertexMap)); 11119566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&bodyEdgeMap)); 11129566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&bodyEdgeGlobalMap)); 11139566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&bodyFaceMap)); 1114c1cad2e7SMatthew G. Knepley 1115c1cad2e7SMatthew G. Knepley *newdm = dm; 1116c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 1117c1cad2e7SMatthew G. Knepley } 1118c1cad2e7SMatthew G. Knepley 1119d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateEGADS_Tess_Internal(MPI_Comm comm, ego context, ego model, DM *newdm) 1120d71ae5a4SJacob Faibussowitsch { 1121c1cad2e7SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel; 1122c1cad2e7SMatthew G. Knepley /* EGADSLite variables */ 1123c1cad2e7SMatthew G. Knepley ego geom, *bodies, *fobjs; 1124c1cad2e7SMatthew G. Knepley int b, oclass, mtype, nbodies, *senses; 1125c1cad2e7SMatthew G. Knepley int totalNumTris = 0, totalNumPoints = 0; 1126c1cad2e7SMatthew G. Knepley double boundBox[6] = {0., 0., 0., 0., 0., 0.}, tessSize; 1127c1cad2e7SMatthew G. Knepley /* PETSc variables */ 1128c1cad2e7SMatthew G. Knepley DM dm; 1129c1cad2e7SMatthew G. Knepley PetscHMapI pointIndexStartMap = NULL, triIndexStartMap = NULL, pTypeLabelMap = NULL, pIndexLabelMap = NULL; 1130c1cad2e7SMatthew G. Knepley PetscHMapI pBodyIndexLabelMap = NULL, triFaceIDLabelMap = NULL, triBodyIDLabelMap = NULL; 1131c1cad2e7SMatthew G. Knepley PetscInt dim = -1, cdim = -1, numCorners = 0, counter = 0; 1132c1cad2e7SMatthew G. Knepley PetscInt *cells = NULL; 1133c1cad2e7SMatthew G. Knepley const PetscInt *cone = NULL; 1134c1cad2e7SMatthew G. Knepley PetscReal *coords = NULL; 1135c1cad2e7SMatthew G. Knepley PetscMPIInt rank; 1136c1cad2e7SMatthew G. Knepley 1137c1cad2e7SMatthew G. Knepley PetscFunctionBeginUser; 11389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1139c5853193SPierre Jolivet if (rank == 0) { 1140c1cad2e7SMatthew G. Knepley // --------------------------------------------------------------------------------------------------- 1141c1cad2e7SMatthew G. Knepley // Generate Petsc Plex from EGADSlite created Tessellation of geometry 1142c1cad2e7SMatthew G. Knepley // --------------------------------------------------------------------------------------------------- 1143c1cad2e7SMatthew G. Knepley 1144*d5b43468SJose E. Roman // Calculate cell and vertex sizes 11459566063dSJacob Faibussowitsch PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses)); 1146c1cad2e7SMatthew G. Knepley 11479566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&pointIndexStartMap)); 11489566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&triIndexStartMap)); 11499566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&pTypeLabelMap)); 11509566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&pIndexLabelMap)); 11519566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&pBodyIndexLabelMap)); 11529566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&triFaceIDLabelMap)); 11539566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&triBodyIDLabelMap)); 1154c1cad2e7SMatthew G. Knepley 1155c1cad2e7SMatthew G. Knepley /* Create Tessellation of Bodies */ 1156c1cad2e7SMatthew G. Knepley ego tessArray[nbodies]; 1157c1cad2e7SMatthew G. Knepley 1158c1cad2e7SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 1159c1cad2e7SMatthew G. Knepley ego body = bodies[b]; 1160c1cad2e7SMatthew G. Knepley double params[3] = {0.0, 0.0, 0.0}; // Parameters for Tessellation 1161c1cad2e7SMatthew G. Knepley int Nf, bodyNumPoints = 0, bodyNumTris = 0; 1162c1cad2e7SMatthew G. Knepley PetscHashIter PISiter, TISiter; 1163c1cad2e7SMatthew G. Knepley PetscBool PISfound, TISfound; 1164c1cad2e7SMatthew G. Knepley 1165c1cad2e7SMatthew G. Knepley /* Store Start Index for each Body's Point and Tris */ 11669566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(pointIndexStartMap, b, &PISiter, &PISfound)); 11679566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(triIndexStartMap, b, &TISiter, &TISfound)); 1168c1cad2e7SMatthew G. Knepley 11699566063dSJacob Faibussowitsch if (!PISfound) PetscCall(PetscHMapISet(pointIndexStartMap, b, totalNumPoints)); 11709566063dSJacob Faibussowitsch if (!TISfound) PetscCall(PetscHMapISet(triIndexStartMap, b, totalNumTris)); 1171c1cad2e7SMatthew G. Knepley 1172c1cad2e7SMatthew G. Knepley /* Calculate Tessellation parameters based on Bounding Box */ 1173c1cad2e7SMatthew G. Knepley /* Get Bounding Box Dimensions of the BODY */ 11749566063dSJacob Faibussowitsch PetscCall(EG_getBoundingBox(body, boundBox)); 1175c1cad2e7SMatthew G. Knepley tessSize = boundBox[3] - boundBox[0]; 1176c1cad2e7SMatthew G. Knepley if (tessSize < boundBox[4] - boundBox[1]) tessSize = boundBox[4] - boundBox[1]; 1177c1cad2e7SMatthew G. Knepley if (tessSize < boundBox[5] - boundBox[2]) tessSize = boundBox[5] - boundBox[2]; 1178c1cad2e7SMatthew G. Knepley 1179c1cad2e7SMatthew G. Knepley // TODO :: May want to give users tessellation parameter options // 1180c1cad2e7SMatthew G. Knepley params[0] = 0.0250 * tessSize; 1181c1cad2e7SMatthew G. Knepley params[1] = 0.0075 * tessSize; 1182c1cad2e7SMatthew G. Knepley params[2] = 15.0; 1183c1cad2e7SMatthew G. Knepley 11849566063dSJacob Faibussowitsch PetscCall(EG_makeTessBody(body, params, &tessArray[b])); 1185c1cad2e7SMatthew G. Knepley 11869566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 1187c1cad2e7SMatthew G. Knepley 1188c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 1189c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 1190c1cad2e7SMatthew G. Knepley int len, fID, ntris; 1191c1cad2e7SMatthew G. Knepley const int *ptype, *pindex, *ptris, *ptric; 1192c1cad2e7SMatthew G. Knepley const double *pxyz, *puv; 1193c1cad2e7SMatthew G. Knepley 1194c1cad2e7SMatthew G. Knepley // Get Face ID // 1195c1cad2e7SMatthew G. Knepley fID = EG_indexBodyTopo(body, face); 1196c1cad2e7SMatthew G. Knepley 1197c1cad2e7SMatthew G. Knepley // Checkout the Surface Tessellation // 11989566063dSJacob Faibussowitsch PetscCall(EG_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric)); 1199c1cad2e7SMatthew G. Knepley 1200c1cad2e7SMatthew G. Knepley // Determine total number of triangle cells in the tessellation // 1201c1cad2e7SMatthew G. Knepley bodyNumTris += (int)ntris; 1202c1cad2e7SMatthew G. Knepley 1203c1cad2e7SMatthew G. Knepley // Check out the point index and coordinate // 1204c1cad2e7SMatthew G. Knepley for (int p = 0; p < len; ++p) { 1205c1cad2e7SMatthew G. Knepley int global; 1206c1cad2e7SMatthew G. Knepley 12079566063dSJacob Faibussowitsch PetscCall(EG_localToGlobal(tessArray[b], fID, p + 1, &global)); 1208c1cad2e7SMatthew G. Knepley 1209c1cad2e7SMatthew G. Knepley // Determine the total number of points in the tessellation // 1210c1cad2e7SMatthew G. Knepley bodyNumPoints = PetscMax(bodyNumPoints, global); 1211c1cad2e7SMatthew G. Knepley } 1212c1cad2e7SMatthew G. Knepley } 1213c1cad2e7SMatthew G. Knepley EG_free(fobjs); 1214c1cad2e7SMatthew G. Knepley 1215c1cad2e7SMatthew G. Knepley totalNumPoints += bodyNumPoints; 1216c1cad2e7SMatthew G. Knepley totalNumTris += bodyNumTris; 1217c1cad2e7SMatthew G. Knepley } 1218c5853193SPierre Jolivet //} - Original End of (rank == 0) 1219c1cad2e7SMatthew G. Knepley 1220c1cad2e7SMatthew G. Knepley dim = 2; 1221c1cad2e7SMatthew G. Knepley cdim = 3; 1222c1cad2e7SMatthew G. Knepley numCorners = 3; 1223c1cad2e7SMatthew G. Knepley //PetscInt counter = 0; 1224c1cad2e7SMatthew G. Knepley 1225c1cad2e7SMatthew G. Knepley /* NEED TO DEFINE MATRICES/VECTORS TO STORE GEOM REFERENCE DATA */ 1226c1cad2e7SMatthew G. Knepley /* Fill in below and use to define DMLabels after DMPlex creation */ 12279566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(totalNumPoints * cdim, &coords, totalNumTris * numCorners, &cells)); 1228c1cad2e7SMatthew G. Knepley 1229c1cad2e7SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 1230c1cad2e7SMatthew G. Knepley ego body = bodies[b]; 1231c1cad2e7SMatthew G. Knepley int Nf; 1232c1cad2e7SMatthew G. Knepley PetscInt pointIndexStart; 1233c1cad2e7SMatthew G. Knepley PetscHashIter PISiter; 1234c1cad2e7SMatthew G. Knepley PetscBool PISfound; 1235c1cad2e7SMatthew G. Knepley 12369566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(pointIndexStartMap, b, &PISiter, &PISfound)); 123728b400f6SJacob Faibussowitsch PetscCheck(PISfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in pointIndexStartMap", b); 12389566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(pointIndexStartMap, b, &pointIndexStart)); 1239c1cad2e7SMatthew G. Knepley 12409566063dSJacob Faibussowitsch PetscCall(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 1241c1cad2e7SMatthew G. Knepley 1242c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 1243c1cad2e7SMatthew G. Knepley /* Get Face Object */ 1244c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 1245c1cad2e7SMatthew G. Knepley int len, fID, ntris; 1246c1cad2e7SMatthew G. Knepley const int *ptype, *pindex, *ptris, *ptric; 1247c1cad2e7SMatthew G. Knepley const double *pxyz, *puv; 1248c1cad2e7SMatthew G. Knepley 1249c1cad2e7SMatthew G. Knepley /* Get Face ID */ 1250c1cad2e7SMatthew G. Knepley fID = EG_indexBodyTopo(body, face); 1251c1cad2e7SMatthew G. Knepley 1252c1cad2e7SMatthew G. Knepley /* Checkout the Surface Tessellation */ 12539566063dSJacob Faibussowitsch PetscCall(EG_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric)); 1254c1cad2e7SMatthew G. Knepley 1255c1cad2e7SMatthew G. Knepley /* Check out the point index and coordinate */ 1256c1cad2e7SMatthew G. Knepley for (int p = 0; p < len; ++p) { 1257c1cad2e7SMatthew G. Knepley int global; 1258c1cad2e7SMatthew G. Knepley PetscHashIter PTLiter, PILiter, PBLiter; 1259c1cad2e7SMatthew G. Knepley PetscBool PTLfound, PILfound, PBLfound; 1260c1cad2e7SMatthew G. Knepley 12619566063dSJacob Faibussowitsch PetscCall(EG_localToGlobal(tessArray[b], fID, p + 1, &global)); 1262c1cad2e7SMatthew G. Knepley 1263c1cad2e7SMatthew G. Knepley /* Set the coordinates array for DAG */ 1264c1cad2e7SMatthew G. Knepley coords[((global - 1 + pointIndexStart) * 3) + 0] = pxyz[(p * 3) + 0]; 1265c1cad2e7SMatthew G. Knepley coords[((global - 1 + pointIndexStart) * 3) + 1] = pxyz[(p * 3) + 1]; 1266c1cad2e7SMatthew G. Knepley coords[((global - 1 + pointIndexStart) * 3) + 2] = pxyz[(p * 3) + 2]; 1267c1cad2e7SMatthew G. Knepley 1268c1cad2e7SMatthew G. Knepley /* Store Geometry Label Information for DMLabel assignment later */ 12699566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(pTypeLabelMap, global - 1 + pointIndexStart, &PTLiter, &PTLfound)); 12709566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(pIndexLabelMap, global - 1 + pointIndexStart, &PILiter, &PILfound)); 12719566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(pBodyIndexLabelMap, global - 1 + pointIndexStart, &PBLiter, &PBLfound)); 1272c1cad2e7SMatthew G. Knepley 12739566063dSJacob Faibussowitsch if (!PTLfound) PetscCall(PetscHMapISet(pTypeLabelMap, global - 1 + pointIndexStart, ptype[p])); 12749566063dSJacob Faibussowitsch if (!PILfound) PetscCall(PetscHMapISet(pIndexLabelMap, global - 1 + pointIndexStart, pindex[p])); 12759566063dSJacob Faibussowitsch if (!PBLfound) PetscCall(PetscHMapISet(pBodyIndexLabelMap, global - 1 + pointIndexStart, b)); 1276c1cad2e7SMatthew G. Knepley 12779566063dSJacob Faibussowitsch if (ptype[p] < 0) PetscCall(PetscHMapISet(pIndexLabelMap, global - 1 + pointIndexStart, fID)); 1278c1cad2e7SMatthew G. Knepley } 1279c1cad2e7SMatthew G. Knepley 1280c1cad2e7SMatthew G. Knepley for (int t = 0; t < (int)ntris; ++t) { 1281c1cad2e7SMatthew G. Knepley int global, globalA, globalB; 1282c1cad2e7SMatthew G. Knepley PetscHashIter TFLiter, TBLiter; 1283c1cad2e7SMatthew G. Knepley PetscBool TFLfound, TBLfound; 1284c1cad2e7SMatthew G. Knepley 12859566063dSJacob Faibussowitsch PetscCall(EG_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 0], &global)); 1286c1cad2e7SMatthew G. Knepley cells[(counter * 3) + 0] = global - 1 + pointIndexStart; 1287c1cad2e7SMatthew G. Knepley 12889566063dSJacob Faibussowitsch PetscCall(EG_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 1], &globalA)); 1289c1cad2e7SMatthew G. Knepley cells[(counter * 3) + 1] = globalA - 1 + pointIndexStart; 1290c1cad2e7SMatthew G. Knepley 12919566063dSJacob Faibussowitsch PetscCall(EG_localToGlobal(tessArray[b], fID, ptris[(t * 3) + 2], &globalB)); 1292c1cad2e7SMatthew G. Knepley cells[(counter * 3) + 2] = globalB - 1 + pointIndexStart; 1293c1cad2e7SMatthew G. Knepley 12949566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(triFaceIDLabelMap, counter, &TFLiter, &TFLfound)); 12959566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(triBodyIDLabelMap, counter, &TBLiter, &TBLfound)); 1296c1cad2e7SMatthew G. Knepley 12979566063dSJacob Faibussowitsch if (!TFLfound) PetscCall(PetscHMapISet(triFaceIDLabelMap, counter, fID)); 12989566063dSJacob Faibussowitsch if (!TBLfound) PetscCall(PetscHMapISet(triBodyIDLabelMap, counter, b)); 1299c1cad2e7SMatthew G. Knepley 1300c1cad2e7SMatthew G. Knepley counter += 1; 1301c1cad2e7SMatthew G. Knepley } 1302c1cad2e7SMatthew G. Knepley } 1303c1cad2e7SMatthew G. Knepley EG_free(fobjs); 1304c1cad2e7SMatthew G. Knepley } 1305c1cad2e7SMatthew G. Knepley } 1306c1cad2e7SMatthew G. Knepley 1307c1cad2e7SMatthew G. Knepley //Build DMPlex 13089566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, totalNumTris, totalNumPoints, numCorners, PETSC_TRUE, cells, cdim, coords, &dm)); 13099566063dSJacob Faibussowitsch PetscCall(PetscFree2(coords, cells)); 1310c1cad2e7SMatthew G. Knepley 1311c1cad2e7SMatthew G. Knepley // Embed EGADS model in DM 1312c1cad2e7SMatthew G. Knepley { 1313c1cad2e7SMatthew G. Knepley PetscContainer modelObj, contextObj; 1314c1cad2e7SMatthew G. Knepley 13159566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &modelObj)); 13169566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(modelObj, model)); 13179566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Model", (PetscObject)modelObj)); 13189566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&modelObj)); 1319c1cad2e7SMatthew G. Knepley 13209566063dSJacob Faibussowitsch PetscCall(PetscContainerCreate(PETSC_COMM_SELF, &contextObj)); 13219566063dSJacob Faibussowitsch PetscCall(PetscContainerSetPointer(contextObj, context)); 13229566063dSJacob Faibussowitsch PetscCall(PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private)); 13239566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "EGADS Context", (PetscObject)contextObj)); 13249566063dSJacob Faibussowitsch PetscCall(PetscContainerDestroy(&contextObj)); 1325c1cad2e7SMatthew G. Knepley } 1326c1cad2e7SMatthew G. Knepley 1327c1cad2e7SMatthew G. Knepley // Label Points 13289566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Body ID")); 13299566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 13309566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Face ID")); 13319566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 13329566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Edge ID")); 13339566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel)); 13349566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "EGADS Vertex ID")); 13359566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel)); 1336c1cad2e7SMatthew G. Knepley 1337c1cad2e7SMatthew G. Knepley /* Get Number of DAG Nodes at each level */ 1338c1cad2e7SMatthew G. Knepley int fStart, fEnd, eStart, eEnd, nStart, nEnd; 1339c1cad2e7SMatthew G. Knepley 13409566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &fStart, &fEnd)); 13419566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 1, &eStart, &eEnd)); 13429566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 2, &nStart, &nEnd)); 1343c1cad2e7SMatthew G. Knepley 1344c1cad2e7SMatthew G. Knepley /* Set DMLabels for NODES */ 1345c1cad2e7SMatthew G. Knepley for (int n = nStart; n < nEnd; ++n) { 1346c1cad2e7SMatthew G. Knepley int pTypeVal, pIndexVal, pBodyVal; 1347c1cad2e7SMatthew G. Knepley PetscHashIter PTLiter, PILiter, PBLiter; 1348c1cad2e7SMatthew G. Knepley PetscBool PTLfound, PILfound, PBLfound; 1349c1cad2e7SMatthew G. Knepley 1350c1cad2e7SMatthew G. Knepley //Converted to Hash Tables 13519566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(pTypeLabelMap, n - nStart, &PTLiter, &PTLfound)); 135228b400f6SJacob Faibussowitsch PetscCheck(PTLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pTypeLabelMap", n); 13539566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(pTypeLabelMap, n - nStart, &pTypeVal)); 1354c1cad2e7SMatthew G. Knepley 13559566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(pIndexLabelMap, n - nStart, &PILiter, &PILfound)); 135628b400f6SJacob Faibussowitsch PetscCheck(PILfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pIndexLabelMap", n); 13579566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(pIndexLabelMap, n - nStart, &pIndexVal)); 1358c1cad2e7SMatthew G. Knepley 13599566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(pBodyIndexLabelMap, n - nStart, &PBLiter, &PBLfound)); 136028b400f6SJacob Faibussowitsch PetscCheck(PBLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pBodyLabelMap", n); 13619566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(pBodyIndexLabelMap, n - nStart, &pBodyVal)); 1362c1cad2e7SMatthew G. Knepley 13639566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, n, pBodyVal)); 13649566063dSJacob Faibussowitsch if (pTypeVal == 0) PetscCall(DMLabelSetValue(vertexLabel, n, pIndexVal)); 13659566063dSJacob Faibussowitsch if (pTypeVal > 0) PetscCall(DMLabelSetValue(edgeLabel, n, pIndexVal)); 13669566063dSJacob Faibussowitsch if (pTypeVal < 0) PetscCall(DMLabelSetValue(faceLabel, n, pIndexVal)); 1367c1cad2e7SMatthew G. Knepley } 1368c1cad2e7SMatthew G. Knepley 1369c1cad2e7SMatthew G. Knepley /* Set DMLabels for Edges - Based on the DMLabels of the EDGE's NODES */ 1370c1cad2e7SMatthew G. Knepley for (int e = eStart; e < eEnd; ++e) { 1371c1cad2e7SMatthew G. Knepley int bodyID_0, vertexID_0, vertexID_1, edgeID_0, edgeID_1, faceID_0, faceID_1; 1372c1cad2e7SMatthew G. Knepley 13739566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, e, &cone)); 13749566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, cone[0], &bodyID_0)); // Do I need to check the other end? 13759566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(vertexLabel, cone[0], &vertexID_0)); 13769566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(vertexLabel, cone[1], &vertexID_1)); 13779566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(edgeLabel, cone[0], &edgeID_0)); 13789566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(edgeLabel, cone[1], &edgeID_1)); 13799566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(faceLabel, cone[0], &faceID_0)); 13809566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(faceLabel, cone[1], &faceID_1)); 1381c1cad2e7SMatthew G. Knepley 13829566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, e, bodyID_0)); 1383c1cad2e7SMatthew G. Knepley 13849566063dSJacob Faibussowitsch if (edgeID_0 == edgeID_1) PetscCall(DMLabelSetValue(edgeLabel, e, edgeID_0)); 13859566063dSJacob Faibussowitsch else if (vertexID_0 > 0 && edgeID_1 > 0) PetscCall(DMLabelSetValue(edgeLabel, e, edgeID_1)); 13869566063dSJacob Faibussowitsch else if (vertexID_1 > 0 && edgeID_0 > 0) PetscCall(DMLabelSetValue(edgeLabel, e, edgeID_0)); 1387c1cad2e7SMatthew G. Knepley else { /* Do Nothing */ } 1388c1cad2e7SMatthew G. Knepley } 1389c1cad2e7SMatthew G. Knepley 1390c1cad2e7SMatthew G. Knepley /* Set DMLabels for Cells */ 1391c1cad2e7SMatthew G. Knepley for (int f = fStart; f < fEnd; ++f) { 1392c1cad2e7SMatthew G. Knepley int edgeID_0; 1393c1cad2e7SMatthew G. Knepley PetscInt triBodyVal, triFaceVal; 1394c1cad2e7SMatthew G. Knepley PetscHashIter TFLiter, TBLiter; 1395c1cad2e7SMatthew G. Knepley PetscBool TFLfound, TBLfound; 1396c1cad2e7SMatthew G. Knepley 1397c1cad2e7SMatthew G. Knepley // Convert to Hash Table 13989566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(triFaceIDLabelMap, f - fStart, &TFLiter, &TFLfound)); 139928b400f6SJacob Faibussowitsch PetscCheck(TFLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in triFaceIDLabelMap", f); 14009566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(triFaceIDLabelMap, f - fStart, &triFaceVal)); 1401c1cad2e7SMatthew G. Knepley 14029566063dSJacob Faibussowitsch PetscCall(PetscHMapIFind(triBodyIDLabelMap, f - fStart, &TBLiter, &TBLfound)); 140328b400f6SJacob Faibussowitsch PetscCheck(TBLfound, PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in triBodyIDLabelMap", f); 14049566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(triBodyIDLabelMap, f - fStart, &triBodyVal)); 1405c1cad2e7SMatthew G. Knepley 14069566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, f, triBodyVal)); 14079566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(faceLabel, f, triFaceVal)); 1408c1cad2e7SMatthew G. Knepley 1409c1cad2e7SMatthew G. Knepley /* Finish Labeling previously unlabeled DMPlex Edges - Assumes Triangular Cell (3 Edges Max) */ 14109566063dSJacob Faibussowitsch PetscCall(DMPlexGetCone(dm, f, &cone)); 1411c1cad2e7SMatthew G. Knepley 1412c1cad2e7SMatthew G. Knepley for (int jj = 0; jj < 3; ++jj) { 14139566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(edgeLabel, cone[jj], &edgeID_0)); 1414c1cad2e7SMatthew G. Knepley 1415c1cad2e7SMatthew G. Knepley if (edgeID_0 < 0) { 14169566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(bodyLabel, cone[jj], triBodyVal)); 14179566063dSJacob Faibussowitsch PetscCall(DMLabelSetValue(faceLabel, cone[jj], triFaceVal)); 1418c1cad2e7SMatthew G. Knepley } 1419c1cad2e7SMatthew G. Knepley } 1420c1cad2e7SMatthew G. Knepley } 1421c1cad2e7SMatthew G. Knepley 1422c1cad2e7SMatthew G. Knepley *newdm = dm; 1423c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 1424c1cad2e7SMatthew G. Knepley } 14257bee2925SMatthew Knepley #endif 14267bee2925SMatthew Knepley 1427c1cad2e7SMatthew G. Knepley /*@ 1428c1cad2e7SMatthew 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. 1429c1cad2e7SMatthew G. Knepley 1430c1cad2e7SMatthew G. Knepley Collective on dm 1431c1cad2e7SMatthew G. Knepley 1432c1cad2e7SMatthew G. Knepley Input Parameter: 1433c1cad2e7SMatthew G. Knepley . dm - The uninflated DM object representing the mesh 1434c1cad2e7SMatthew G. Knepley 1435c1cad2e7SMatthew G. Knepley Output Parameter: 1436c1cad2e7SMatthew G. Knepley . dm - The inflated DM object representing the mesh 1437c1cad2e7SMatthew G. Knepley 1438c1cad2e7SMatthew G. Knepley Level: intermediate 1439c1cad2e7SMatthew G. Knepley 1440db781477SPatrick Sanan .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateEGADS()` 1441c1cad2e7SMatthew G. Knepley @*/ 1442d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexInflateToGeomModel(DM dm) 1443d71ae5a4SJacob Faibussowitsch { 1444c1cad2e7SMatthew G. Knepley #if defined(PETSC_HAVE_EGADS) 1445c1cad2e7SMatthew G. Knepley /* EGADS Variables */ 1446c1cad2e7SMatthew G. Knepley ego model, geom, body, face, edge; 1447c1cad2e7SMatthew G. Knepley ego *bodies; 1448c1cad2e7SMatthew G. Knepley int Nb, oclass, mtype, *senses; 1449c1cad2e7SMatthew G. Knepley double result[3]; 1450c1cad2e7SMatthew G. Knepley /* PETSc Variables */ 1451c1cad2e7SMatthew G. Knepley DM cdm; 1452c1cad2e7SMatthew G. Knepley PetscContainer modelObj; 1453c1cad2e7SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel; 1454c1cad2e7SMatthew G. Knepley Vec coordinates; 1455c1cad2e7SMatthew G. Knepley PetscScalar *coords; 1456c1cad2e7SMatthew G. Knepley PetscInt bodyID, faceID, edgeID, vertexID; 1457c1cad2e7SMatthew G. Knepley PetscInt cdim, d, vStart, vEnd, v; 1458c1cad2e7SMatthew G. Knepley #endif 1459c1cad2e7SMatthew G. Knepley 1460c1cad2e7SMatthew G. Knepley PetscFunctionBegin; 1461c1cad2e7SMatthew G. Knepley #if defined(PETSC_HAVE_EGADS) 14629566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "EGADS Model", (PetscObject *)&modelObj)); 1463c1cad2e7SMatthew G. Knepley if (!modelObj) PetscFunctionReturn(0); 14649566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &cdim)); 14659566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 14669566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 14679566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 14689566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 14699566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel)); 14709566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel)); 1471c1cad2e7SMatthew G. Knepley 14729566063dSJacob Faibussowitsch PetscCall(PetscContainerGetPointer(modelObj, (void **)&model)); 14739566063dSJacob Faibussowitsch PetscCall(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 1474c1cad2e7SMatthew G. Knepley 14759566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 14769566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(coordinates, &coords)); 1477c1cad2e7SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 1478c1cad2e7SMatthew G. Knepley PetscScalar *vcoords; 1479c1cad2e7SMatthew G. Knepley 14809566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(bodyLabel, v, &bodyID)); 14819566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(faceLabel, v, &faceID)); 14829566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(edgeLabel, v, &edgeID)); 14839566063dSJacob Faibussowitsch PetscCall(DMLabelGetValue(vertexLabel, v, &vertexID)); 1484c1cad2e7SMatthew G. Knepley 148563a3b9bcSJacob Faibussowitsch PetscCheck(bodyID < Nb, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %" PetscInt_FMT " is not in [0, %d)", bodyID, Nb); 1486c1cad2e7SMatthew G. Knepley body = bodies[bodyID]; 1487c1cad2e7SMatthew G. Knepley 14889566063dSJacob Faibussowitsch PetscCall(DMPlexPointLocalRef(cdm, v, coords, (void *)&vcoords)); 1489c1cad2e7SMatthew G. Knepley if (edgeID > 0) { 1490c1cad2e7SMatthew G. Knepley /* Snap to EDGE at nearest location */ 1491c1cad2e7SMatthew G. Knepley double params[1]; 14929566063dSJacob Faibussowitsch PetscCall(EG_objectBodyTopo(body, EDGE, edgeID, &edge)); 14939566063dSJacob Faibussowitsch PetscCall(EG_invEvaluate(edge, vcoords, params, result)); // Get (x,y,z) of nearest point on EDGE 1494c1cad2e7SMatthew G. Knepley for (d = 0; d < cdim; ++d) vcoords[d] = result[d]; 1495c1cad2e7SMatthew G. Knepley } else if (faceID > 0) { 1496c1cad2e7SMatthew G. Knepley /* Snap to FACE at nearest location */ 1497c1cad2e7SMatthew G. Knepley double params[2]; 14989566063dSJacob Faibussowitsch PetscCall(EG_objectBodyTopo(body, FACE, faceID, &face)); 14999566063dSJacob Faibussowitsch PetscCall(EG_invEvaluate(face, vcoords, params, result)); // Get (x,y,z) of nearest point on FACE 1500c1cad2e7SMatthew G. Knepley for (d = 0; d < cdim; ++d) vcoords[d] = result[d]; 1501c1cad2e7SMatthew G. Knepley } 1502c1cad2e7SMatthew G. Knepley } 15039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(coordinates, &coords)); 1504c1cad2e7SMatthew G. Knepley /* Clear out global coordinates */ 15056858538eSMatthew G. Knepley PetscCall(VecDestroy(&dm->coordinates[0].x)); 1506c1cad2e7SMatthew G. Knepley #endif 1507c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 1508c1cad2e7SMatthew G. Knepley } 1509c1cad2e7SMatthew G. Knepley 15107bee2925SMatthew Knepley /*@C 1511c1cad2e7SMatthew G. Knepley DMPlexCreateEGADSFromFile - Create a DMPlex mesh from an EGADS, IGES, or STEP file. 15127bee2925SMatthew Knepley 15137bee2925SMatthew Knepley Collective 15147bee2925SMatthew Knepley 15157bee2925SMatthew Knepley Input Parameters: 15167bee2925SMatthew Knepley + comm - The MPI communicator 1517c1cad2e7SMatthew G. Knepley - filename - The name of the EGADS, IGES, or STEP file 15187bee2925SMatthew Knepley 15197bee2925SMatthew Knepley Output Parameter: 15207bee2925SMatthew Knepley . dm - The DM object representing the mesh 15217bee2925SMatthew Knepley 15227bee2925SMatthew Knepley Level: beginner 15237bee2925SMatthew Knepley 1524db781477SPatrick Sanan .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateEGADS()`, `DMPlexCreateEGADSLiteFromFile()` 15257bee2925SMatthew Knepley @*/ 1526d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateEGADSFromFile(MPI_Comm comm, const char filename[], DM *dm) 1527d71ae5a4SJacob Faibussowitsch { 15287bee2925SMatthew Knepley PetscMPIInt rank; 15297bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS) 15307bee2925SMatthew Knepley ego context = NULL, model = NULL; 15317bee2925SMatthew Knepley #endif 1532c1cad2e7SMatthew G. Knepley PetscBool printModel = PETSC_FALSE, tessModel = PETSC_FALSE, newModel = PETSC_FALSE; 15337bee2925SMatthew Knepley 15347bee2925SMatthew Knepley PetscFunctionBegin; 15357bee2925SMatthew Knepley PetscValidCharPointer(filename, 2); 15369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_print_model", &printModel, NULL)); 15379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_tess_model", &tessModel, NULL)); 15389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_new_model", &newModel, NULL)); 15399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 15407bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS) 1541dd400576SPatrick Sanan if (rank == 0) { 15429566063dSJacob Faibussowitsch PetscCall(EG_open(&context)); 15439566063dSJacob Faibussowitsch PetscCall(EG_loadModel(context, 0, filename, &model)); 15449566063dSJacob Faibussowitsch if (printModel) PetscCall(DMPlexEGADSPrintModel_Internal(model)); 15457bee2925SMatthew Knepley } 15469566063dSJacob Faibussowitsch if (tessModel) PetscCall(DMPlexCreateEGADS_Tess_Internal(comm, context, model, dm)); 15479566063dSJacob Faibussowitsch else if (newModel) PetscCall(DMPlexCreateEGADS(comm, context, model, dm)); 15489566063dSJacob Faibussowitsch else PetscCall(DMPlexCreateEGADS_Internal(comm, context, model, dm)); 1549c03d1177SJose E. Roman PetscFunctionReturn(0); 15507bee2925SMatthew Knepley #else 1551c1cad2e7SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "This method requires EGADS support. Reconfigure using --download-egads"); 15527bee2925SMatthew Knepley #endif 15537bee2925SMatthew Knepley } 1554