1a8ededdfSMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 27bee2925SMatthew Knepley #include <petsc/private/hashmapi.h> 3a8ededdfSMatthew G. Knepley 4a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS 5a8ededdfSMatthew G. Knepley #include <egads.h> 6a8ededdfSMatthew G. Knepley #endif 7a8ededdfSMatthew G. Knepley 8a8ededdfSMatthew G. Knepley /* We need to understand how to natively parse STEP files. There seems to be only one open source implementation of 9a8ededdfSMatthew G. Knepley the STEP parser contained in the OpenCASCADE package. It is enough to make a strong man weep: 10a8ededdfSMatthew G. Knepley 11a8ededdfSMatthew G. Knepley https://github.com/tpaviot/oce/tree/master/src/STEPControl 12a8ededdfSMatthew G. Knepley 13a8ededdfSMatthew G. Knepley The STEP, and inner EXPRESS, formats are ISO standards, so they are documented 14a8ededdfSMatthew G. Knepley 15a8ededdfSMatthew G. Knepley https://stackoverflow.com/questions/26774037/documentation-or-specification-for-step-and-stp-files 16a8ededdfSMatthew G. Knepley http://stepmod.sourceforge.net/express_model_spec/ 17a8ededdfSMatthew G. Knepley 18a8ededdfSMatthew G. Knepley but again it seems that there has been a deliberate effort at obfuscation, probably to raise the bar for entrants. 19a8ededdfSMatthew G. Knepley */ 20a8ededdfSMatthew G. Knepley 21c1cad2e7SMatthew G. Knepley #ifdef PETSC_HAVE_EGADS 22c1cad2e7SMatthew G. Knepley PETSC_INTERN PetscErrorCode DMPlexSnapToGeomModel_EGADS_Internal(DM, PetscInt, ego, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]); 23c1cad2e7SMatthew G. Knepley PETSC_INTERN PetscErrorCode DMPlexSnapToGeomModel_EGADSLite_Internal(DM, PetscInt, ego, PetscInt, PetscInt, PetscInt, const PetscScalar[], PetscScalar[]); 24c1cad2e7SMatthew G. Knepley 25c1cad2e7SMatthew G. Knepley PetscErrorCode DMPlexSnapToGeomModel_EGADS_Internal(DM dm, PetscInt p, ego model, PetscInt bodyID, PetscInt faceID, PetscInt edgeID, const PetscScalar mcoords[], PetscScalar gcoords[]) 26c1cad2e7SMatthew G. Knepley { 27c1cad2e7SMatthew G. Knepley DM cdm; 28c1cad2e7SMatthew G. Knepley ego *bodies; 29c1cad2e7SMatthew G. Knepley ego geom, body, obj; 30a5b23f4aSJose E. Roman /* result has to hold derivatives, along with the value */ 31c1cad2e7SMatthew G. Knepley double params[3], result[18], paramsV[16*3], resultV[16*3], range[4]; 32c1cad2e7SMatthew G. Knepley int Nb, oclass, mtype, *senses, peri; 33c1cad2e7SMatthew G. Knepley Vec coordinatesLocal; 34c1cad2e7SMatthew G. Knepley PetscScalar *coords = NULL; 35c1cad2e7SMatthew G. Knepley PetscInt Nv, v, Np = 0, pm; 36c1cad2e7SMatthew G. Knepley PetscInt dE, d; 37c1cad2e7SMatthew G. Knepley 38c1cad2e7SMatthew G. Knepley PetscFunctionBeginHot; 395f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dm, &cdm)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDim(dm, &dE)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordinatesLocal)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 432c71b3e2SJacob Faibussowitsch PetscCheckFalse(bodyID >= Nb,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %D is not in [0, %d)", bodyID, Nb); 44c1cad2e7SMatthew G. Knepley body = bodies[bodyID]; 45c1cad2e7SMatthew G. Knepley 465f80ce2aSJacob Faibussowitsch if (edgeID >= 0) {CHKERRQ(EG_objectBodyTopo(body, EDGE, edgeID, &obj)); Np = 1;} 475f80ce2aSJacob Faibussowitsch else if (faceID >= 0) {CHKERRQ(EG_objectBodyTopo(body, FACE, faceID, &obj)); Np = 2;} 48c1cad2e7SMatthew G. Knepley else { 49c1cad2e7SMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 50c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 51c1cad2e7SMatthew G. Knepley } 52c1cad2e7SMatthew G. Knepley 53c1cad2e7SMatthew G. Knepley /* Calculate parameters (t or u,v) for vertices */ 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecGetClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords)); 55c1cad2e7SMatthew G. Knepley Nv /= dE; 56c1cad2e7SMatthew G. Knepley if (Nv == 1) { 575f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords)); 58c1cad2e7SMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 59c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 60c1cad2e7SMatthew G. Knepley } 612c71b3e2SJacob Faibussowitsch PetscCheckFalse(Nv > 16,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot handle %D coordinates associated to point %D", Nv, p); 62c1cad2e7SMatthew G. Knepley 63c1cad2e7SMatthew G. Knepley /* Correct EGADSlite 2pi bug when calculating nearest point on Periodic Surfaces */ 645f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getRange(obj, range, &peri)); 65c1cad2e7SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 665f80ce2aSJacob Faibussowitsch CHKERRQ(EG_invEvaluate(obj, &coords[v*dE], ¶msV[v*3], &resultV[v*3])); 67c1cad2e7SMatthew G. Knepley #if 1 68c1cad2e7SMatthew G. Knepley if (peri > 0) { 69c1cad2e7SMatthew G. Knepley if (paramsV[v*3+0] + 1.e-4 < range[0]) {paramsV[v*3+0] += 2. * PETSC_PI;} 70c1cad2e7SMatthew G. Knepley else if (paramsV[v*3+0] - 1.e-4 > range[1]) {paramsV[v*3+0] -= 2. * PETSC_PI;} 71c1cad2e7SMatthew G. Knepley } 72c1cad2e7SMatthew G. Knepley if (peri > 1) { 73c1cad2e7SMatthew G. Knepley if (paramsV[v*3+1] + 1.e-4 < range[2]) {paramsV[v*3+1] += 2. * PETSC_PI;} 74c1cad2e7SMatthew G. Knepley else if (paramsV[v*3+1] - 1.e-4 > range[3]) {paramsV[v*3+1] -= 2. * PETSC_PI;} 75c1cad2e7SMatthew G. Knepley } 76c1cad2e7SMatthew G. Knepley #endif 77c1cad2e7SMatthew G. Knepley } 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexVecRestoreClosure(cdm, NULL, coordinatesLocal, p, &Nv, &coords)); 79c1cad2e7SMatthew G. Knepley /* Calculate parameters (t or u,v) for new vertex at edge midpoint */ 80c1cad2e7SMatthew G. Knepley for (pm = 0; pm < Np; ++pm) { 81c1cad2e7SMatthew G. Knepley params[pm] = 0.; 82c1cad2e7SMatthew G. Knepley for (v = 0; v < Nv; ++v) params[pm] += paramsV[v*3+pm]; 83c1cad2e7SMatthew G. Knepley params[pm] /= Nv; 84c1cad2e7SMatthew G. Knepley } 852c71b3e2SJacob Faibussowitsch PetscCheckFalse((params[0] < range[0]) || (params[0] > range[1]),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D had bad interpolation", p); 862c71b3e2SJacob Faibussowitsch PetscCheckFalse(Np > 1 && ((params[1] < range[2]) || (params[1] > range[3])),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D had bad interpolation", p); 87c1cad2e7SMatthew G. Knepley /* Put coordinates for new vertex in result[] */ 885f80ce2aSJacob Faibussowitsch CHKERRQ(EG_evaluate(obj, params, result)); 89c1cad2e7SMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = result[d]; 90c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 91c1cad2e7SMatthew G. Knepley } 92c1cad2e7SMatthew G. Knepley #endif 93c1cad2e7SMatthew G. Knepley 94a8ededdfSMatthew G. Knepley /*@ 95a8ededdfSMatthew 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. 96a8ededdfSMatthew G. Knepley 97a8ededdfSMatthew G. Knepley Not collective 98a8ededdfSMatthew G. Knepley 99a8ededdfSMatthew G. Knepley Input Parameters: 100a8ededdfSMatthew G. Knepley + dm - The DMPlex object 101a8ededdfSMatthew G. Knepley . p - The mesh point 102d410b0cfSMatthew G. Knepley . dE - The coordinate dimension 103a8ededdfSMatthew G. Knepley - mcoords - A coordinate point lying on the mesh point 104a8ededdfSMatthew G. Knepley 105a8ededdfSMatthew G. Knepley Output Parameter: 106a8ededdfSMatthew G. Knepley . gcoords - The closest coordinate point on the geometry model associated with 'p' to the given point 107a8ededdfSMatthew G. Knepley 108d410b0cfSMatthew 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. 109a8ededdfSMatthew G. Knepley 110a8ededdfSMatthew G. Knepley Level: intermediate 111a8ededdfSMatthew G. Knepley 112a8ededdfSMatthew G. Knepley .seealso: DMRefine(), DMPlexCreate(), DMPlexSetRefinementUniform() 113a8ededdfSMatthew G. Knepley @*/ 114d410b0cfSMatthew G. Knepley PetscErrorCode DMPlexSnapToGeomModel(DM dm, PetscInt p, PetscInt dE, const PetscScalar mcoords[], PetscScalar gcoords[]) 115a8ededdfSMatthew G. Knepley { 116d410b0cfSMatthew G. Knepley PetscInt d; 117a8ededdfSMatthew G. Knepley 118c1cad2e7SMatthew G. Knepley PetscFunctionBeginHot; 119a8ededdfSMatthew G. Knepley #ifdef PETSC_HAVE_EGADS 120c1cad2e7SMatthew G. Knepley { 121c1cad2e7SMatthew G. Knepley DM_Plex *plex = (DM_Plex *) dm->data; 122c1cad2e7SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel; 123c1cad2e7SMatthew G. Knepley PetscInt bodyID, faceID, edgeID; 124c1cad2e7SMatthew G. Knepley PetscContainer modelObj; 125c1cad2e7SMatthew G. Knepley ego model; 126c1cad2e7SMatthew G. Knepley PetscBool islite = PETSC_FALSE; 127c1cad2e7SMatthew G. Knepley 1285f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel)); 131c1cad2e7SMatthew G. Knepley if (!bodyLabel || !faceLabel || !edgeLabel || plex->ignoreModel) { 132f0fcf0adSMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 133a8ededdfSMatthew G. Knepley PetscFunctionReturn(0); 134a8ededdfSMatthew G. Knepley } 1355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj)); 136c1cad2e7SMatthew G. Knepley if (!modelObj) { 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectQuery((PetscObject) dm, "EGADSLite Model", (PetscObject *) &modelObj)); 138c1cad2e7SMatthew G. Knepley islite = PETSC_TRUE; 139c1cad2e7SMatthew G. Knepley } 140c1cad2e7SMatthew G. Knepley if (!modelObj) { 141c1cad2e7SMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 142c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 143c1cad2e7SMatthew G. Knepley } 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerGetPointer(modelObj, (void **) &model)); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(bodyLabel, p, &bodyID)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(faceLabel, p, &faceID)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(edgeLabel, p, &edgeID)); 148c1cad2e7SMatthew G. Knepley /* Allows for "Connective" Plex Edges present in models with multiple non-touching Entities */ 149c1cad2e7SMatthew G. Knepley if (bodyID < 0) { 150f0fcf0adSMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 151f0fcf0adSMatthew G. Knepley PetscFunctionReturn(0); 152a8ededdfSMatthew G. Knepley } 1535f80ce2aSJacob Faibussowitsch if (islite) CHKERRQ(DMPlexSnapToGeomModel_EGADSLite_Internal(dm, p, model, bodyID, faceID, edgeID, mcoords, gcoords)); 1545f80ce2aSJacob Faibussowitsch else CHKERRQ(DMPlexSnapToGeomModel_EGADS_Internal(dm, p, model, bodyID, faceID, edgeID, mcoords, gcoords)); 155f0fcf0adSMatthew G. Knepley } 156a8ededdfSMatthew G. Knepley #else 157f0fcf0adSMatthew G. Knepley for (d = 0; d < dE; ++d) gcoords[d] = mcoords[d]; 158a8ededdfSMatthew G. Knepley #endif 159a8ededdfSMatthew G. Knepley PetscFunctionReturn(0); 160a8ededdfSMatthew G. Knepley } 1617bee2925SMatthew Knepley 1627bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS) 163c1cad2e7SMatthew G. Knepley static PetscErrorCode DMPlexEGADSPrintModel_Internal(ego model) 1647bee2925SMatthew Knepley { 1657bee2925SMatthew Knepley ego geom, *bodies, *objs, *nobjs, *mobjs, *lobjs; 1667bee2925SMatthew Knepley int oclass, mtype, *senses; 1677bee2925SMatthew Knepley int Nb, b; 1687bee2925SMatthew Knepley 1697bee2925SMatthew Knepley PetscFunctionBeginUser; 1707bee2925SMatthew Knepley /* test bodyTopo functions */ 1715f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 1725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Number of BODIES (nbodies): %d \n", Nb)); 1737bee2925SMatthew Knepley 1747bee2925SMatthew Knepley for (b = 0; b < Nb; ++b) { 1757bee2925SMatthew Knepley ego body = bodies[b]; 1767bee2925SMatthew Knepley int id, Nsh, Nf, Nl, l, Ne, e, Nv, v; 1777bee2925SMatthew Knepley 1787bee2925SMatthew Knepley /* Output Basic Model Topology */ 1795f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, SHELL, &Nsh, &objs)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Number of SHELLS: %d \n", Nsh)); 1817bee2925SMatthew Knepley EG_free(objs); 1827bee2925SMatthew Knepley 1835f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, FACE, &Nf, &objs)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Number of FACES: %d \n", Nf)); 1857bee2925SMatthew Knepley EG_free(objs); 1867bee2925SMatthew Knepley 1875f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs)); 1885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Number of LOOPS: %d \n", Nl)); 1897bee2925SMatthew Knepley 1905f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, EDGE, &Ne, &objs)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Number of EDGES: %d \n", Ne)); 1927bee2925SMatthew Knepley EG_free(objs); 1937bee2925SMatthew Knepley 1945f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, NODE, &Nv, &objs)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Number of NODES: %d \n", Nv)); 1967bee2925SMatthew Knepley EG_free(objs); 1977bee2925SMatthew Knepley 1987bee2925SMatthew Knepley for (l = 0; l < Nl; ++l) { 1997bee2925SMatthew Knepley ego loop = lobjs[l]; 2007bee2925SMatthew Knepley 2017bee2925SMatthew Knepley id = EG_indexBodyTopo(body, loop); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " LOOP ID: %d\n", id)); 2037bee2925SMatthew Knepley 2047bee2925SMatthew Knepley /* Get EDGE info which associated with the current LOOP */ 2055f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses)); 2067bee2925SMatthew Knepley 2077bee2925SMatthew Knepley for (e = 0; e < Ne; ++e) { 2087bee2925SMatthew Knepley ego edge = objs[e]; 2097bee2925SMatthew Knepley double range[4] = {0., 0., 0., 0.}; 2107bee2925SMatthew Knepley double point[3] = {0., 0., 0.}; 211266cfabeSMatthew G. Knepley double params[3] = {0., 0., 0.}; 212266cfabeSMatthew G. Knepley double result[18]; 2137bee2925SMatthew Knepley int peri; 2147bee2925SMatthew Knepley 2155f80ce2aSJacob Faibussowitsch id = EG_indexBodyTopo(body, edge); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " EDGE ID: %d (%d)\n", id, e)); 2177bee2925SMatthew Knepley 2185f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getRange(edge, range, &peri)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " Range = %lf, %lf, %lf, %lf \n", range[0], range[1], range[2], range[3])); 2207bee2925SMatthew Knepley 2217bee2925SMatthew Knepley /* Get NODE info which associated with the current EDGE */ 2225f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 223266cfabeSMatthew G. Knepley if (mtype == DEGENERATE) { 2245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " EDGE %d is DEGENERATE \n", id)); 225266cfabeSMatthew G. Knepley } else { 226266cfabeSMatthew G. Knepley params[0] = range[0]; 2275f80ce2aSJacob Faibussowitsch CHKERRQ(EG_evaluate(edge, params, result)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " between (%lf, %lf, %lf)", result[0], result[1], result[2])); 229266cfabeSMatthew G. Knepley params[0] = range[1]; 2305f80ce2aSJacob Faibussowitsch CHKERRQ(EG_evaluate(edge, params, result)); 2315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " and (%lf, %lf, %lf)\n", result[0], result[1], result[2])); 232266cfabeSMatthew G. Knepley } 2337bee2925SMatthew Knepley 2347bee2925SMatthew Knepley for (v = 0; v < Nv; ++v) { 2357bee2925SMatthew Knepley ego vertex = nobjs[v]; 2367bee2925SMatthew Knepley double limits[4]; 2377bee2925SMatthew Knepley int dummy; 2387bee2925SMatthew Knepley 2395f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses)); 2407bee2925SMatthew Knepley id = EG_indexBodyTopo(body, vertex); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " NODE ID: %d \n", id)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " (x, y, z) = (%lf, %lf, %lf) \n", limits[0], limits[1], limits[2])); 2437bee2925SMatthew Knepley 2447bee2925SMatthew Knepley point[0] = point[0] + limits[0]; 2457bee2925SMatthew Knepley point[1] = point[1] + limits[1]; 2467bee2925SMatthew Knepley point[2] = point[2] + limits[2]; 2477bee2925SMatthew Knepley } 2487bee2925SMatthew Knepley } 2497bee2925SMatthew Knepley } 2507bee2925SMatthew Knepley EG_free(lobjs); 2517bee2925SMatthew Knepley } 2527bee2925SMatthew Knepley PetscFunctionReturn(0); 2537bee2925SMatthew Knepley } 2547bee2925SMatthew Knepley 2557bee2925SMatthew Knepley static PetscErrorCode DMPlexEGADSDestroy_Private(void *context) 2567bee2925SMatthew Knepley { 2577bee2925SMatthew Knepley if (context) EG_close((ego) context); 2587bee2925SMatthew Knepley return 0; 2597bee2925SMatthew Knepley } 2607bee2925SMatthew Knepley 261c1cad2e7SMatthew G. Knepley static PetscErrorCode DMPlexCreateEGADS_Internal(MPI_Comm comm, ego context, ego model, DM *newdm) 2627bee2925SMatthew Knepley { 263c1cad2e7SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel; 2647bee2925SMatthew Knepley PetscInt cStart, cEnd, c; 2657bee2925SMatthew Knepley /* EGADSLite variables */ 2667bee2925SMatthew Knepley ego geom, *bodies, *objs, *nobjs, *mobjs, *lobjs; 2677bee2925SMatthew Knepley int oclass, mtype, nbodies, *senses; 2687bee2925SMatthew Knepley int b; 2697bee2925SMatthew Knepley /* PETSc variables */ 2707bee2925SMatthew Knepley DM dm; 2717bee2925SMatthew Knepley PetscHMapI edgeMap = NULL; 2727bee2925SMatthew 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; 2737bee2925SMatthew Knepley PetscInt *cells = NULL, *cone = NULL; 2747bee2925SMatthew Knepley PetscReal *coords = NULL; 2757bee2925SMatthew Knepley PetscMPIInt rank; 2767bee2925SMatthew Knepley 2777bee2925SMatthew Knepley PetscFunctionBegin; 2785f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 279dd400576SPatrick Sanan if (rank == 0) { 280266cfabeSMatthew G. Knepley const PetscInt debug = 0; 281266cfabeSMatthew G. Knepley 2827bee2925SMatthew Knepley /* --------------------------------------------------------------------------------------------------- 2837bee2925SMatthew Knepley Generate Petsc Plex 2847bee2925SMatthew Knepley Get all Nodes in model, record coordinates in a correctly formatted array 2857bee2925SMatthew Knepley Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array 2867bee2925SMatthew Knepley We need to uniformly refine the initial geometry to guarantee a valid mesh 2877bee2925SMatthew Knepley */ 2887bee2925SMatthew Knepley 2897bee2925SMatthew Knepley /* Calculate cell and vertex sizes */ 2905f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses)); 2915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&edgeMap)); 2927bee2925SMatthew Knepley numEdges = 0; 2937bee2925SMatthew Knepley for (b = 0; b < nbodies; ++b) { 2947bee2925SMatthew Knepley ego body = bodies[b]; 2957bee2925SMatthew Knepley int id, Nl, l, Nv, v; 2967bee2925SMatthew Knepley 2975f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs)); 2987bee2925SMatthew Knepley for (l = 0; l < Nl; ++l) { 2997bee2925SMatthew Knepley ego loop = lobjs[l]; 3007bee2925SMatthew Knepley int Ner = 0, Ne, e, Nc; 3017bee2925SMatthew Knepley 3025f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses)); 3037bee2925SMatthew Knepley for (e = 0; e < Ne; ++e) { 3047bee2925SMatthew Knepley ego edge = objs[e]; 3057bee2925SMatthew Knepley int Nv, id; 3067bee2925SMatthew Knepley PetscHashIter iter; 3077bee2925SMatthew Knepley PetscBool found; 3087bee2925SMatthew Knepley 3095f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 3107bee2925SMatthew Knepley if (mtype == DEGENERATE) continue; 3115f80ce2aSJacob Faibussowitsch id = EG_indexBodyTopo(body, edge); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(edgeMap, id-1, &iter, &found)); 3135f80ce2aSJacob Faibussowitsch if (!found) CHKERRQ(PetscHMapISet(edgeMap, id-1, numEdges++)); 3147bee2925SMatthew Knepley ++Ner; 3157bee2925SMatthew Knepley } 3167bee2925SMatthew Knepley if (Ner == 2) {Nc = 2;} 3177bee2925SMatthew Knepley else if (Ner == 3) {Nc = 4;} 3187bee2925SMatthew Knepley else if (Ner == 4) {Nc = 8; ++numQuads;} 31998921bdaSJacob Faibussowitsch else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot support loop with %d edges", Ner); 3207bee2925SMatthew Knepley numCells += Nc; 3217bee2925SMatthew Knepley newCells += Nc-1; 3227bee2925SMatthew Knepley maxCorners = PetscMax(Ner*2+1, maxCorners); 3237bee2925SMatthew Knepley } 3245f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs)); 3257bee2925SMatthew Knepley for (v = 0; v < Nv; ++v) { 3267bee2925SMatthew Knepley ego vertex = nobjs[v]; 3277bee2925SMatthew Knepley 3287bee2925SMatthew Knepley id = EG_indexBodyTopo(body, vertex); 3297bee2925SMatthew Knepley /* TODO: Instead of assuming contiguous ids, we could use a hash table */ 3307bee2925SMatthew Knepley numVertices = PetscMax(id, numVertices); 3317bee2925SMatthew Knepley } 3327bee2925SMatthew Knepley EG_free(lobjs); 3337bee2925SMatthew Knepley EG_free(nobjs); 3347bee2925SMatthew Knepley } 3355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGetSize(edgeMap, &numEdges)); 3367bee2925SMatthew Knepley newVertices = numEdges + numQuads; 3377bee2925SMatthew Knepley numVertices += newVertices; 3387bee2925SMatthew Knepley 3397bee2925SMatthew Knepley dim = 2; /* Assume 3D Models :: Need to update to handle 2D Models in the future */ 3407bee2925SMatthew Knepley cdim = 3; /* Assume 3D Models :: Need to update to handle 2D Models in the future */ 3417bee2925SMatthew Knepley numCorners = 3; /* Split cells into triangles */ 3425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(numVertices*cdim, &coords, numCells*numCorners, &cells, maxCorners, &cone)); 3437bee2925SMatthew Knepley 3447bee2925SMatthew Knepley /* Get vertex coordinates */ 3457bee2925SMatthew Knepley for (b = 0; b < nbodies; ++b) { 3467bee2925SMatthew Knepley ego body = bodies[b]; 3477bee2925SMatthew Knepley int id, Nv, v; 3487bee2925SMatthew Knepley 3495f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs)); 3507bee2925SMatthew Knepley for (v = 0; v < Nv; ++v) { 3517bee2925SMatthew Knepley ego vertex = nobjs[v]; 3527bee2925SMatthew Knepley double limits[4]; 3537bee2925SMatthew Knepley int dummy; 3547bee2925SMatthew Knepley 3555f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses)); 3565f80ce2aSJacob Faibussowitsch id = EG_indexBodyTopo(body, vertex); 3577bee2925SMatthew Knepley coords[(id-1)*cdim+0] = limits[0]; 3587bee2925SMatthew Knepley coords[(id-1)*cdim+1] = limits[1]; 3597bee2925SMatthew Knepley coords[(id-1)*cdim+2] = limits[2]; 3607bee2925SMatthew Knepley } 3617bee2925SMatthew Knepley EG_free(nobjs); 3627bee2925SMatthew Knepley } 3635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIClear(edgeMap)); 3647bee2925SMatthew Knepley fOff = numVertices - newVertices + numEdges; 3657bee2925SMatthew Knepley numEdges = 0; 3667bee2925SMatthew Knepley numQuads = 0; 3677bee2925SMatthew Knepley for (b = 0; b < nbodies; ++b) { 3687bee2925SMatthew Knepley ego body = bodies[b]; 3697bee2925SMatthew Knepley int Nl, l; 3707bee2925SMatthew Knepley 3715f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs)); 3727bee2925SMatthew Knepley for (l = 0; l < Nl; ++l) { 3737bee2925SMatthew Knepley ego loop = lobjs[l]; 3747bee2925SMatthew Knepley int lid, Ner = 0, Ne, e; 3757bee2925SMatthew Knepley 3765f80ce2aSJacob Faibussowitsch lid = EG_indexBodyTopo(body, loop); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses)); 3787bee2925SMatthew Knepley for (e = 0; e < Ne; ++e) { 3797bee2925SMatthew Knepley ego edge = objs[e]; 3807bee2925SMatthew Knepley int eid, Nv; 3817bee2925SMatthew Knepley PetscHashIter iter; 3827bee2925SMatthew Knepley PetscBool found; 3837bee2925SMatthew Knepley 3845f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 3857bee2925SMatthew Knepley if (mtype == DEGENERATE) continue; 3867bee2925SMatthew Knepley ++Ner; 3875f80ce2aSJacob Faibussowitsch eid = EG_indexBodyTopo(body, edge); 3885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(edgeMap, eid-1, &iter, &found)); 3897bee2925SMatthew Knepley if (!found) { 3907bee2925SMatthew Knepley PetscInt v = numVertices - newVertices + numEdges; 3917bee2925SMatthew Knepley double range[4], params[3] = {0., 0., 0.}, result[18]; 3927bee2925SMatthew Knepley int periodic[2]; 3937bee2925SMatthew Knepley 3945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapISet(edgeMap, eid-1, numEdges++)); 3955f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getRange(edge, range, periodic)); 3967bee2925SMatthew Knepley params[0] = 0.5*(range[0] + range[1]); 3975f80ce2aSJacob Faibussowitsch CHKERRQ(EG_evaluate(edge, params, result)); 3987bee2925SMatthew Knepley coords[v*cdim+0] = result[0]; 3997bee2925SMatthew Knepley coords[v*cdim+1] = result[1]; 4007bee2925SMatthew Knepley coords[v*cdim+2] = result[2]; 4017bee2925SMatthew Knepley } 4027bee2925SMatthew Knepley } 4037bee2925SMatthew Knepley if (Ner == 4) { 4047bee2925SMatthew Knepley PetscInt v = fOff + numQuads++; 405266cfabeSMatthew G. Knepley ego *fobjs, face; 4067bee2925SMatthew Knepley double range[4], params[3] = {0., 0., 0.}, result[18]; 407266cfabeSMatthew G. Knepley int Nf, fid, periodic[2]; 4087bee2925SMatthew Knepley 4095f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs)); 410266cfabeSMatthew G. Knepley face = fobjs[0]; 4115f80ce2aSJacob Faibussowitsch fid = EG_indexBodyTopo(body, face); 4122c71b3e2SJacob Faibussowitsch PetscCheckFalse(Nf != 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Loop %d has %d faces, instead of 1 (%d)", lid-1, Nf, fid); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getRange(face, range, periodic)); 4147bee2925SMatthew Knepley params[0] = 0.5*(range[0] + range[1]); 4157bee2925SMatthew Knepley params[1] = 0.5*(range[2] + range[3]); 4165f80ce2aSJacob Faibussowitsch CHKERRQ(EG_evaluate(face, params, result)); 4177bee2925SMatthew Knepley coords[v*cdim+0] = result[0]; 4187bee2925SMatthew Knepley coords[v*cdim+1] = result[1]; 4197bee2925SMatthew Knepley coords[v*cdim+2] = result[2]; 4207bee2925SMatthew Knepley } 4217bee2925SMatthew Knepley } 4227bee2925SMatthew Knepley } 4232c71b3e2SJacob Faibussowitsch PetscCheckFalse(numEdges + numQuads != newVertices,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of new vertices %D != %D previous count", numEdges + numQuads, newVertices); 4247bee2925SMatthew Knepley 4257bee2925SMatthew Knepley /* Get cell vertices by traversing loops */ 4267bee2925SMatthew Knepley numQuads = 0; 4277bee2925SMatthew Knepley cOff = 0; 4287bee2925SMatthew Knepley for (b = 0; b < nbodies; ++b) { 4297bee2925SMatthew Knepley ego body = bodies[b]; 4307bee2925SMatthew Knepley int id, Nl, l; 4317bee2925SMatthew Knepley 4325f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs)); 4337bee2925SMatthew Knepley for (l = 0; l < Nl; ++l) { 4347bee2925SMatthew Knepley ego loop = lobjs[l]; 4357bee2925SMatthew Knepley int lid, Ner = 0, Ne, e, nc = 0, c, Nt, t; 4367bee2925SMatthew Knepley 4375f80ce2aSJacob Faibussowitsch lid = EG_indexBodyTopo(body, loop); 4385f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses)); 4397bee2925SMatthew Knepley 4407bee2925SMatthew Knepley for (e = 0; e < Ne; ++e) { 4417bee2925SMatthew Knepley ego edge = objs[e]; 4427bee2925SMatthew Knepley int points[3]; 4437bee2925SMatthew Knepley int eid, Nv, v, tmp; 4447bee2925SMatthew Knepley 4457bee2925SMatthew Knepley eid = EG_indexBodyTopo(body, edge); 4465f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 447266cfabeSMatthew G. Knepley if (mtype == DEGENERATE) continue; 448266cfabeSMatthew G. Knepley else ++Ner; 4492c71b3e2SJacob Faibussowitsch PetscCheckFalse(Nv != 2,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Edge %d has %d vertices != 2", eid, Nv); 4507bee2925SMatthew Knepley 4517bee2925SMatthew Knepley for (v = 0; v < Nv; ++v) { 4527bee2925SMatthew Knepley ego vertex = nobjs[v]; 4537bee2925SMatthew Knepley 4547bee2925SMatthew Knepley id = EG_indexBodyTopo(body, vertex); 4557bee2925SMatthew Knepley points[v*2] = id-1; 4567bee2925SMatthew Knepley } 4577bee2925SMatthew Knepley { 4587bee2925SMatthew Knepley PetscInt edgeNum; 4597bee2925SMatthew Knepley 4605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(edgeMap, eid-1, &edgeNum)); 4617bee2925SMatthew Knepley points[1] = numVertices - newVertices + edgeNum; 4627bee2925SMatthew Knepley } 4637bee2925SMatthew Knepley /* EGADS loops are not oriented, but seem to be in order, so we must piece them together */ 4647bee2925SMatthew Knepley if (!nc) { 4657bee2925SMatthew Knepley for (v = 0; v < Nv+1; ++v) cone[nc++] = points[v]; 4667bee2925SMatthew Knepley } else { 4677bee2925SMatthew Knepley if (cone[nc-1] == points[0]) {cone[nc++] = points[1]; if (cone[0] != points[2]) cone[nc++] = points[2];} 4687bee2925SMatthew Knepley else if (cone[nc-1] == points[2]) {cone[nc++] = points[1]; if (cone[0] != points[0]) cone[nc++] = points[0];} 4697bee2925SMatthew Knepley else if (cone[nc-3] == points[0]) {tmp = cone[nc-3]; cone[nc-3] = cone[nc-1]; cone[nc-1] = tmp; cone[nc++] = points[1]; if (cone[0] != points[2]) cone[nc++] = points[2];} 4707bee2925SMatthew Knepley else if (cone[nc-3] == points[2]) {tmp = cone[nc-3]; cone[nc-3] = cone[nc-1]; cone[nc-1] = tmp; cone[nc++] = points[1]; if (cone[0] != points[0]) cone[nc++] = points[0];} 47198921bdaSJacob Faibussowitsch else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %d does not match its predecessor", eid); 4727bee2925SMatthew Knepley } 4737bee2925SMatthew Knepley } 4742c71b3e2SJacob Faibussowitsch PetscCheckFalse(nc != 2*Ner,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %D != %D", nc, 2*Ner); 4757bee2925SMatthew Knepley if (Ner == 4) {cone[nc++] = numVertices - newVertices + numEdges + numQuads++;} 4762c71b3e2SJacob Faibussowitsch PetscCheckFalse(nc > maxCorners,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %D > %D max", nc, maxCorners); 4777bee2925SMatthew Knepley /* Triangulate the loop */ 4787bee2925SMatthew Knepley switch (Ner) { 4797bee2925SMatthew Knepley case 2: /* Bi-Segment -> 2 triangles */ 4807bee2925SMatthew Knepley Nt = 2; 4817bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[0]; 4827bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[1]; 4837bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[2]; 4847bee2925SMatthew Knepley ++cOff; 4857bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[0]; 4867bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[2]; 4877bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[3]; 4887bee2925SMatthew Knepley ++cOff; 4897bee2925SMatthew Knepley break; 4907bee2925SMatthew Knepley case 3: /* Triangle -> 4 triangles */ 4917bee2925SMatthew Knepley Nt = 4; 4927bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[0]; 4937bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[1]; 4947bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[5]; 4957bee2925SMatthew Knepley ++cOff; 4967bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[1]; 4977bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[2]; 4987bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[3]; 4997bee2925SMatthew Knepley ++cOff; 5007bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[5]; 5017bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[3]; 5027bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[4]; 5037bee2925SMatthew Knepley ++cOff; 5047bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[1]; 5057bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[3]; 5067bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[5]; 5077bee2925SMatthew Knepley ++cOff; 5087bee2925SMatthew Knepley break; 5097bee2925SMatthew Knepley case 4: /* Quad -> 8 triangles */ 5107bee2925SMatthew Knepley Nt = 8; 5117bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[0]; 5127bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[1]; 5137bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[7]; 5147bee2925SMatthew Knepley ++cOff; 5157bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[1]; 5167bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[2]; 5177bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[3]; 5187bee2925SMatthew Knepley ++cOff; 5197bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[3]; 5207bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[4]; 5217bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[5]; 5227bee2925SMatthew Knepley ++cOff; 5237bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[5]; 5247bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[6]; 5257bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[7]; 5267bee2925SMatthew Knepley ++cOff; 5277bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[8]; 5287bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[1]; 5297bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[3]; 5307bee2925SMatthew Knepley ++cOff; 5317bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[8]; 5327bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[3]; 5337bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[5]; 5347bee2925SMatthew Knepley ++cOff; 5357bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[8]; 5367bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[5]; 5377bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[7]; 5387bee2925SMatthew Knepley ++cOff; 5397bee2925SMatthew Knepley cells[cOff*numCorners+0] = cone[8]; 5407bee2925SMatthew Knepley cells[cOff*numCorners+1] = cone[7]; 5417bee2925SMatthew Knepley cells[cOff*numCorners+2] = cone[1]; 5427bee2925SMatthew Knepley ++cOff; 5437bee2925SMatthew Knepley break; 54498921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d edges, which we do not support", lid, Ner); 5457bee2925SMatthew Knepley } 546266cfabeSMatthew G. Knepley if (debug) { 5477bee2925SMatthew Knepley for (t = 0; t < Nt; ++t) { 5485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, " LOOP Corner NODEs Triangle %D (", t)); 5497bee2925SMatthew Knepley for (c = 0; c < numCorners; ++c) { 5505f80ce2aSJacob Faibussowitsch if (c > 0) CHKERRQ(PetscPrintf(PETSC_COMM_SELF, ", ")); 5515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, "%D", cells[(cOff-Nt+t)*numCorners+c])); 5527bee2925SMatthew Knepley } 5535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF, ")\n")); 5547bee2925SMatthew Knepley } 5557bee2925SMatthew Knepley } 556266cfabeSMatthew G. Knepley } 5577bee2925SMatthew Knepley EG_free(lobjs); 5587bee2925SMatthew Knepley } 5597bee2925SMatthew Knepley } 5602c71b3e2SJacob Faibussowitsch PetscCheckFalse(cOff != numCells,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count of total cells %D != %D previous count", cOff, numCells); 5615f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numVertices, numCorners, PETSC_TRUE, cells, cdim, coords, &dm)); 5625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(coords, cells, cone)); 5635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(dm, " Total Number of Unique Cells = %D (%D)\n", numCells, newCells)); 5645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(dm, " Total Number of Unique Vertices = %D (%D)\n", numVertices, newVertices)); 5657bee2925SMatthew Knepley /* Embed EGADS model in DM */ 5667bee2925SMatthew Knepley { 5677bee2925SMatthew Knepley PetscContainer modelObj, contextObj; 5687bee2925SMatthew Knepley 5695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerCreate(PETSC_COMM_SELF, &modelObj)); 5705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetPointer(modelObj, model)); 5715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose((PetscObject) dm, "EGADS Model", (PetscObject) modelObj)); 5725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerDestroy(&modelObj)); 5737bee2925SMatthew Knepley 5745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerCreate(PETSC_COMM_SELF, &contextObj)); 5755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetPointer(contextObj, context)); 5765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private)); 5775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose((PetscObject) dm, "EGADS Context", (PetscObject) contextObj)); 5785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerDestroy(&contextObj)); 5797bee2925SMatthew Knepley } 5807bee2925SMatthew Knepley /* Label points */ 5815f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Body ID")); 5825f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 5835f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Face ID")); 5845f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 5855f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Edge ID")); 5865f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel)); 5875f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Vertex ID")); 5885f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel)); 5897bee2925SMatthew Knepley cOff = 0; 5907bee2925SMatthew Knepley for (b = 0; b < nbodies; ++b) { 5917bee2925SMatthew Knepley ego body = bodies[b]; 5927bee2925SMatthew Knepley int id, Nl, l; 5937bee2925SMatthew Knepley 5945f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs)); 5957bee2925SMatthew Knepley for (l = 0; l < Nl; ++l) { 5967bee2925SMatthew Knepley ego loop = lobjs[l]; 5977bee2925SMatthew Knepley ego *fobjs; 5987bee2925SMatthew Knepley int lid, Nf, fid, Ner = 0, Ne, e, Nt = 0, t; 5997bee2925SMatthew Knepley 6005f80ce2aSJacob Faibussowitsch lid = EG_indexBodyTopo(body, loop); 6015f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs)); 6022c71b3e2SJacob Faibussowitsch PetscCheckFalse(Nf > 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf); 6035f80ce2aSJacob Faibussowitsch fid = EG_indexBodyTopo(body, fobjs[0]); 6047bee2925SMatthew Knepley EG_free(fobjs); 6055f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses)); 6067bee2925SMatthew Knepley for (e = 0; e < Ne; ++e) { 6077bee2925SMatthew Knepley ego edge = objs[e]; 6087bee2925SMatthew Knepley int eid, Nv, v; 6097bee2925SMatthew Knepley PetscInt points[3], support[2], numEdges, edgeNum; 6107bee2925SMatthew Knepley const PetscInt *edges; 6117bee2925SMatthew Knepley 6127bee2925SMatthew Knepley eid = EG_indexBodyTopo(body, edge); 6135f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 6147bee2925SMatthew Knepley if (mtype == DEGENERATE) continue; 6157bee2925SMatthew Knepley else ++Ner; 6167bee2925SMatthew Knepley for (v = 0; v < Nv; ++v) { 6177bee2925SMatthew Knepley ego vertex = nobjs[v]; 6187bee2925SMatthew Knepley 6197bee2925SMatthew Knepley id = EG_indexBodyTopo(body, vertex); 6205f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(edgeLabel, numCells + id-1, eid)); 6217bee2925SMatthew Knepley points[v*2] = numCells + id-1; 6227bee2925SMatthew Knepley } 6235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(edgeMap, eid-1, &edgeNum)); 6247bee2925SMatthew Knepley points[1] = numCells + numVertices - newVertices + edgeNum; 6257bee2925SMatthew Knepley 6265f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(edgeLabel, points[1], eid)); 6277bee2925SMatthew Knepley support[0] = points[0]; 6287bee2925SMatthew Knepley support[1] = points[1]; 6295f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetJoin(dm, 2, support, &numEdges, &edges)); 6302c71b3e2SJacob Faibussowitsch PetscCheckFalse(numEdges != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vertices (%D, %D) should only bound 1 edge, not %D", support[0], support[1], numEdges); 6315f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(edgeLabel, edges[0], eid)); 6325f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges)); 6337bee2925SMatthew Knepley support[0] = points[1]; 6347bee2925SMatthew Knepley support[1] = points[2]; 6355f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetJoin(dm, 2, support, &numEdges, &edges)); 6362c71b3e2SJacob Faibussowitsch PetscCheckFalse(numEdges != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vertices (%D, %D) should only bound 1 edge, not %D", support[0], support[1], numEdges); 6375f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(edgeLabel, edges[0], eid)); 6385f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges)); 6397bee2925SMatthew Knepley } 6407bee2925SMatthew Knepley switch (Ner) { 6417bee2925SMatthew Knepley case 2: Nt = 2;break; 6427bee2925SMatthew Knepley case 3: Nt = 4;break; 6437bee2925SMatthew Knepley case 4: Nt = 8;break; 64498921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Loop with %d edges is unsupported", Ner); 6457bee2925SMatthew Knepley } 6467bee2925SMatthew Knepley for (t = 0; t < Nt; ++t) { 6475f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, cOff+t, b)); 6485f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(faceLabel, cOff+t, fid)); 6497bee2925SMatthew Knepley } 6507bee2925SMatthew Knepley cOff += Nt; 6517bee2925SMatthew Knepley } 6527bee2925SMatthew Knepley EG_free(lobjs); 6537bee2925SMatthew Knepley } 6545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIDestroy(&edgeMap)); 6555f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 6567bee2925SMatthew Knepley for (c = cStart; c < cEnd; ++c) { 6577bee2925SMatthew Knepley PetscInt *closure = NULL; 6587bee2925SMatthew Knepley PetscInt clSize, cl, bval, fval; 6597bee2925SMatthew Knepley 6605f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure)); 6615f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(bodyLabel, c, &bval)); 6625f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(faceLabel, c, &fval)); 6637bee2925SMatthew Knepley for (cl = 0; cl < clSize*2; cl += 2) { 6645f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, closure[cl], bval)); 6655f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(faceLabel, closure[cl], fval)); 6667bee2925SMatthew Knepley } 6675f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure)); 6687bee2925SMatthew Knepley } 6697bee2925SMatthew Knepley *newdm = dm; 6707bee2925SMatthew Knepley PetscFunctionReturn(0); 6717bee2925SMatthew Knepley } 672c1cad2e7SMatthew G. Knepley 673c1cad2e7SMatthew G. Knepley static PetscErrorCode DMPlexCreateEGADS(MPI_Comm comm, ego context, ego model, DM *newdm) 674c1cad2e7SMatthew G. Knepley { 675c1cad2e7SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel; 676c1cad2e7SMatthew G. Knepley // EGADS/EGADSLite variables 677c1cad2e7SMatthew G. Knepley ego geom, *bodies, *mobjs, *fobjs, *lobjs, *eobjs, *nobjs; 678c1cad2e7SMatthew G. Knepley ego topRef, prev, next; 679c1cad2e7SMatthew G. Knepley int oclass, mtype, nbodies, *senses, *lSenses, *eSenses; 680c1cad2e7SMatthew G. Knepley int b; 681c1cad2e7SMatthew G. Knepley // PETSc variables 682c1cad2e7SMatthew G. Knepley DM dm; 683c1cad2e7SMatthew G. Knepley PetscHMapI edgeMap = NULL, bodyIndexMap = NULL, bodyVertexMap = NULL, bodyEdgeMap = NULL, bodyFaceMap = NULL, bodyEdgeGlobalMap = NULL; 684c1cad2e7SMatthew G. Knepley PetscInt dim = -1, cdim = -1, numCorners = 0, numVertices = 0, numEdges = 0, numFaces = 0, numCells = 0, edgeCntr = 0; 685c1cad2e7SMatthew G. Knepley PetscInt cellCntr = 0, numPoints = 0; 686c1cad2e7SMatthew G. Knepley PetscInt *cells = NULL; 687c1cad2e7SMatthew G. Knepley const PetscInt *cone = NULL; 688c1cad2e7SMatthew G. Knepley PetscReal *coords = NULL; 689c1cad2e7SMatthew G. Knepley PetscMPIInt rank; 690c1cad2e7SMatthew G. Knepley 691c1cad2e7SMatthew G. Knepley PetscFunctionBeginUser; 6925f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 693c1cad2e7SMatthew G. Knepley if (!rank) { 694c1cad2e7SMatthew G. Knepley // --------------------------------------------------------------------------------------------------- 695c1cad2e7SMatthew G. Knepley // Generate Petsc Plex 696c1cad2e7SMatthew G. Knepley // Get all Nodes in model, record coordinates in a correctly formatted array 697c1cad2e7SMatthew G. Knepley // Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array 698c1cad2e7SMatthew G. Knepley // We need to uniformly refine the initial geometry to guarantee a valid mesh 699c1cad2e7SMatthew G. Knepley 700c1cad2e7SMatthew G. Knepley // Caluculate cell and vertex sizes 7015f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses)); 702c1cad2e7SMatthew G. Knepley 7035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&edgeMap)); 7045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&bodyIndexMap)); 7055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&bodyVertexMap)); 7065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&bodyEdgeMap)); 7075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&bodyEdgeGlobalMap)); 7085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&bodyFaceMap)); 709c1cad2e7SMatthew G. Knepley 710c1cad2e7SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 711c1cad2e7SMatthew G. Knepley ego body = bodies[b]; 712c1cad2e7SMatthew G. Knepley int Nf, Ne, Nv; 713c1cad2e7SMatthew G. Knepley PetscHashIter BIiter, BViter, BEiter, BEGiter, BFiter, EMiter; 714c1cad2e7SMatthew G. Knepley PetscBool BIfound, BVfound, BEfound, BEGfound, BFfound, EMfound; 715c1cad2e7SMatthew G. Knepley 7165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyIndexMap, b, &BIiter, &BIfound)); 7175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound)); 7185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound)); 7195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound)); 7205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound)); 721c1cad2e7SMatthew G. Knepley 7225f80ce2aSJacob Faibussowitsch if (!BIfound) CHKERRQ(PetscHMapISet(bodyIndexMap, b, numFaces + numEdges + numVertices)); 7235f80ce2aSJacob Faibussowitsch if (!BVfound) CHKERRQ(PetscHMapISet(bodyVertexMap, b, numVertices)); 7245f80ce2aSJacob Faibussowitsch if (!BEfound) CHKERRQ(PetscHMapISet(bodyEdgeMap, b, numEdges)); 7255f80ce2aSJacob Faibussowitsch if (!BEGfound) CHKERRQ(PetscHMapISet(bodyEdgeGlobalMap, b, edgeCntr)); 7265f80ce2aSJacob Faibussowitsch if (!BFfound) CHKERRQ(PetscHMapISet(bodyFaceMap, b, numFaces)); 727c1cad2e7SMatthew G. Knepley 7285f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 7295f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs)); 7305f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs)); 731c1cad2e7SMatthew G. Knepley EG_free(fobjs); 732c1cad2e7SMatthew G. Knepley EG_free(eobjs); 733c1cad2e7SMatthew G. Knepley EG_free(nobjs); 734c1cad2e7SMatthew G. Knepley 735c1cad2e7SMatthew G. Knepley // Remove DEGENERATE EDGES from Edge count 7365f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs)); 737c1cad2e7SMatthew G. Knepley int Netemp = 0; 738c1cad2e7SMatthew G. Knepley for (int e = 0; e < Ne; ++e) { 739c1cad2e7SMatthew G. Knepley ego edge = eobjs[e]; 740c1cad2e7SMatthew G. Knepley int eid; 741c1cad2e7SMatthew G. Knepley 7425f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next)); 7435f80ce2aSJacob Faibussowitsch eid = EG_indexBodyTopo(body, edge); 744c1cad2e7SMatthew G. Knepley 7455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(edgeMap, edgeCntr + eid - 1, &EMiter, &EMfound)); 746c1cad2e7SMatthew G. Knepley if (mtype == DEGENERATE) { 7475f80ce2aSJacob Faibussowitsch if (!EMfound) CHKERRQ(PetscHMapISet(edgeMap, edgeCntr + eid - 1, -1)); 748c1cad2e7SMatthew G. Knepley } 749c1cad2e7SMatthew G. Knepley else { 750c1cad2e7SMatthew G. Knepley ++Netemp; 7515f80ce2aSJacob Faibussowitsch if (!EMfound) CHKERRQ(PetscHMapISet(edgeMap, edgeCntr + eid - 1, Netemp)); 752c1cad2e7SMatthew G. Knepley } 753c1cad2e7SMatthew G. Knepley } 754c1cad2e7SMatthew G. Knepley EG_free(eobjs); 755c1cad2e7SMatthew G. Knepley 756c1cad2e7SMatthew G. Knepley // Determine Number of Cells 7575f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 758c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 759c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 760c1cad2e7SMatthew G. Knepley int edgeTemp = 0; 761c1cad2e7SMatthew G. Knepley 7625f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, face, EDGE, &Ne, &eobjs)); 763c1cad2e7SMatthew G. Knepley for (int e = 0; e < Ne; ++e) { 764c1cad2e7SMatthew G. Knepley ego edge = eobjs[e]; 765c1cad2e7SMatthew G. Knepley 7665f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next)); 767c1cad2e7SMatthew G. Knepley if (mtype != DEGENERATE) {++edgeTemp;} 768c1cad2e7SMatthew G. Knepley } 769c1cad2e7SMatthew G. Knepley numCells += (2 * edgeTemp); 770c1cad2e7SMatthew G. Knepley EG_free(eobjs); 771c1cad2e7SMatthew G. Knepley } 772c1cad2e7SMatthew G. Knepley EG_free(fobjs); 773c1cad2e7SMatthew G. Knepley 774c1cad2e7SMatthew G. Knepley numFaces += Nf; 775c1cad2e7SMatthew G. Knepley numEdges += Netemp; 776c1cad2e7SMatthew G. Knepley numVertices += Nv; 777c1cad2e7SMatthew G. Knepley edgeCntr += Ne; 778c1cad2e7SMatthew G. Knepley } 779c1cad2e7SMatthew G. Knepley 780c1cad2e7SMatthew G. Knepley // Set up basic DMPlex parameters 781c1cad2e7SMatthew G. Knepley dim = 2; // Assumes 3D Models :: Need to handle 2D modles in the future 782c1cad2e7SMatthew G. Knepley cdim = 3; // Assumes 3D Models :: Need to update to handle 2D modles in future 783c1cad2e7SMatthew G. Knepley numCorners = 3; // Split Faces into triangles 784c1cad2e7SMatthew G. Knepley numPoints = numVertices + numEdges + numFaces; // total number of coordinate points 785c1cad2e7SMatthew G. Knepley 7865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(numPoints*cdim, &coords, numCells*numCorners, &cells)); 787c1cad2e7SMatthew G. Knepley 788c1cad2e7SMatthew G. Knepley // Get Vertex Coordinates and Set up Cells 789c1cad2e7SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 790c1cad2e7SMatthew G. Knepley ego body = bodies[b]; 791c1cad2e7SMatthew G. Knepley int Nf, Ne, Nv; 792c1cad2e7SMatthew G. Knepley PetscInt bodyVertexIndexStart, bodyEdgeIndexStart, bodyEdgeGlobalIndexStart, bodyFaceIndexStart; 793c1cad2e7SMatthew G. Knepley PetscHashIter BViter, BEiter, BEGiter, BFiter, EMiter; 794c1cad2e7SMatthew G. Knepley PetscBool BVfound, BEfound, BEGfound, BFfound, EMfound; 795c1cad2e7SMatthew G. Knepley 796c1cad2e7SMatthew G. Knepley // Vertices on Current Body 7975f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs)); 798c1cad2e7SMatthew G. Knepley 7995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound)); 800*28b400f6SJacob Faibussowitsch PetscCheck(BVfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyVertexMap", b); 8015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(bodyVertexMap, b, &bodyVertexIndexStart)); 802c1cad2e7SMatthew G. Knepley 803c1cad2e7SMatthew G. Knepley for (int v = 0; v < Nv; ++v) { 804c1cad2e7SMatthew G. Knepley ego vertex = nobjs[v]; 805c1cad2e7SMatthew G. Knepley double limits[4]; 806c1cad2e7SMatthew G. Knepley int id, dummy; 807c1cad2e7SMatthew G. Knepley 8085f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses)); 8095f80ce2aSJacob Faibussowitsch id = EG_indexBodyTopo(body, vertex); 810c1cad2e7SMatthew G. Knepley 811c1cad2e7SMatthew G. Knepley coords[(bodyVertexIndexStart + id - 1)*cdim + 0] = limits[0]; 812c1cad2e7SMatthew G. Knepley coords[(bodyVertexIndexStart + id - 1)*cdim + 1] = limits[1]; 813c1cad2e7SMatthew G. Knepley coords[(bodyVertexIndexStart + id - 1)*cdim + 2] = limits[2]; 814c1cad2e7SMatthew G. Knepley } 815c1cad2e7SMatthew G. Knepley EG_free(nobjs); 816c1cad2e7SMatthew G. Knepley 817c1cad2e7SMatthew G. Knepley // Edge Midpoint Vertices on Current Body 8185f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, EDGE, &Ne, &eobjs)); 819c1cad2e7SMatthew G. Knepley 8205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound)); 821*28b400f6SJacob Faibussowitsch PetscCheck(BEfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeMap", b); 8225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(bodyEdgeMap, b, &bodyEdgeIndexStart)); 823c1cad2e7SMatthew G. Knepley 8245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound)); 825*28b400f6SJacob Faibussowitsch PetscCheck(BEGfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeGlobalMap", b); 8265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(bodyEdgeGlobalMap, b, &bodyEdgeGlobalIndexStart)); 827c1cad2e7SMatthew G. Knepley 828c1cad2e7SMatthew G. Knepley for (int e = 0; e < Ne; ++e) { 829c1cad2e7SMatthew G. Knepley ego edge = eobjs[e]; 830c1cad2e7SMatthew G. Knepley double range[2], avgt[1], cntrPnt[9]; 831c1cad2e7SMatthew G. Knepley int eid, eOffset; 832c1cad2e7SMatthew G. Knepley int periodic; 833c1cad2e7SMatthew G. Knepley 8345f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next)); 835c1cad2e7SMatthew G. Knepley if (mtype == DEGENERATE) {continue;} 836c1cad2e7SMatthew G. Knepley 8375f80ce2aSJacob Faibussowitsch eid = EG_indexBodyTopo(body, edge); 838c1cad2e7SMatthew G. Knepley 839c1cad2e7SMatthew G. Knepley // get relative offset from globalEdgeID Vector 8405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound)); 841*28b400f6SJacob Faibussowitsch PetscCheck(EMfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "Edge %d not found in edgeMap", bodyEdgeGlobalIndexStart + eid - 1); 8425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset)); 843c1cad2e7SMatthew G. Knepley 8445f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getRange(edge, range, &periodic)); 845c1cad2e7SMatthew G. Knepley avgt[0] = (range[0] + range[1]) / 2.; 846c1cad2e7SMatthew G. Knepley 8475f80ce2aSJacob Faibussowitsch CHKERRQ(EG_evaluate(edge, avgt, cntrPnt)); 848c1cad2e7SMatthew G. Knepley coords[(numVertices + bodyEdgeIndexStart + eOffset - 1)*cdim + 0] = cntrPnt[0]; 849c1cad2e7SMatthew G. Knepley coords[(numVertices + bodyEdgeIndexStart + eOffset - 1)*cdim + 1] = cntrPnt[1]; 850c1cad2e7SMatthew G. Knepley coords[(numVertices + bodyEdgeIndexStart + eOffset - 1)*cdim + 2] = cntrPnt[2]; 851c1cad2e7SMatthew G. Knepley } 852c1cad2e7SMatthew G. Knepley EG_free(eobjs); 853c1cad2e7SMatthew G. Knepley 854c1cad2e7SMatthew G. Knepley // Face Midpoint Vertices on Current Body 8555f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 856c1cad2e7SMatthew G. Knepley 8575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound)); 858*28b400f6SJacob Faibussowitsch PetscCheck(BFfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyFaceMap", b); 8595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(bodyFaceMap, b, &bodyFaceIndexStart)); 860c1cad2e7SMatthew G. Knepley 861c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 862c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 863c1cad2e7SMatthew G. Knepley double range[4], avgUV[2], cntrPnt[18]; 864c1cad2e7SMatthew G. Knepley int peri, id; 865c1cad2e7SMatthew G. Knepley 866c1cad2e7SMatthew G. Knepley id = EG_indexBodyTopo(body, face); 8675f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getRange(face, range, &peri)); 868c1cad2e7SMatthew G. Knepley 869c1cad2e7SMatthew G. Knepley avgUV[0] = (range[0] + range[1]) / 2.; 870c1cad2e7SMatthew G. Knepley avgUV[1] = (range[2] + range[3]) / 2.; 8715f80ce2aSJacob Faibussowitsch CHKERRQ(EG_evaluate(face, avgUV, cntrPnt)); 872c1cad2e7SMatthew G. Knepley 873c1cad2e7SMatthew G. Knepley coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1)*cdim + 0] = cntrPnt[0]; 874c1cad2e7SMatthew G. Knepley coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1)*cdim + 1] = cntrPnt[1]; 875c1cad2e7SMatthew G. Knepley coords[(numVertices + numEdges + bodyFaceIndexStart + id - 1)*cdim + 2] = cntrPnt[2]; 876c1cad2e7SMatthew G. Knepley } 877c1cad2e7SMatthew G. Knepley EG_free(fobjs); 878c1cad2e7SMatthew G. Knepley 879c1cad2e7SMatthew G. Knepley // Define Cells :: Note - This could be incorporated in the Face Midpoint Vertices Loop but was kept separate for clarity 8805f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 881c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 882c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 883c1cad2e7SMatthew G. Knepley int fID, midFaceID, midPntID, startID, endID, Nl; 884c1cad2e7SMatthew G. Knepley 8855f80ce2aSJacob Faibussowitsch fID = EG_indexBodyTopo(body, face); 886c1cad2e7SMatthew G. Knepley midFaceID = numVertices + numEdges + bodyFaceIndexStart + fID - 1; 887c1cad2e7SMatthew G. Knepley // Must Traverse Loop to ensure we have all necessary information like the sense (+/- 1) of the edges. 888c1cad2e7SMatthew G. Knepley // TODO :: Only handles single loop faces (No holes). The choices for handling multiloop faces are: 889c1cad2e7SMatthew G. Knepley // 1) Use the DMPlexCreateEGADSFromFile() with the -dm_plex_egads_with_tess = 1 option. 890c1cad2e7SMatthew G. Knepley // This will use a default EGADS tessellation as an initial surface mesh. 891c1cad2e7SMatthew G. Knepley // 2) Create the initial surface mesh via a 2D mesher :: Currently not availble (?future?) 892c1cad2e7SMatthew G. Knepley // May I suggest the XXXX as a starting point? 893c1cad2e7SMatthew G. Knepley 8945f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(face, &geom, &oclass, &mtype, NULL, &Nl, &lobjs, &lSenses)); 895c1cad2e7SMatthew G. Knepley 8962c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 897c1cad2e7SMatthew G. Knepley for (int l = 0; l < Nl; ++l) { 898c1cad2e7SMatthew G. Knepley ego loop = lobjs[l]; 899c1cad2e7SMatthew G. Knepley 9005f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses)); 901c1cad2e7SMatthew G. Knepley for (int e = 0; e < Ne; ++e) { 902c1cad2e7SMatthew G. Knepley ego edge = eobjs[e]; 903c1cad2e7SMatthew G. Knepley int eid, eOffset; 904c1cad2e7SMatthew G. Knepley 9055f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next)); 906c1cad2e7SMatthew G. Knepley eid = EG_indexBodyTopo(body, edge); 907c1cad2e7SMatthew G. Knepley if (mtype == DEGENERATE) { continue; } 908c1cad2e7SMatthew G. Knepley 909c1cad2e7SMatthew G. Knepley // get relative offset from globalEdgeID Vector 9105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound)); 911*28b400f6SJacob 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); 9125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset)); 913c1cad2e7SMatthew G. Knepley 914c1cad2e7SMatthew G. Knepley midPntID = numVertices + bodyEdgeIndexStart + eOffset - 1; 915c1cad2e7SMatthew G. Knepley 9165f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses)); 917c1cad2e7SMatthew G. Knepley 918c1cad2e7SMatthew G. Knepley if (eSenses[e] > 0) { startID = EG_indexBodyTopo(body, nobjs[0]); endID = EG_indexBodyTopo(body, nobjs[1]); } 919c1cad2e7SMatthew G. Knepley else { startID = EG_indexBodyTopo(body, nobjs[1]); endID = EG_indexBodyTopo(body, nobjs[0]); } 920c1cad2e7SMatthew G. Knepley 921c1cad2e7SMatthew G. Knepley // Define 2 Cells per Edge with correct orientation 922c1cad2e7SMatthew G. Knepley cells[cellCntr*numCorners + 0] = midFaceID; 923c1cad2e7SMatthew G. Knepley cells[cellCntr*numCorners + 1] = bodyVertexIndexStart + startID - 1; 924c1cad2e7SMatthew G. Knepley cells[cellCntr*numCorners + 2] = midPntID; 925c1cad2e7SMatthew G. Knepley 926c1cad2e7SMatthew G. Knepley cells[cellCntr*numCorners + 3] = midFaceID; 927c1cad2e7SMatthew G. Knepley cells[cellCntr*numCorners + 4] = midPntID; 928c1cad2e7SMatthew G. Knepley cells[cellCntr*numCorners + 5] = bodyVertexIndexStart + endID - 1; 929c1cad2e7SMatthew G. Knepley 930c1cad2e7SMatthew G. Knepley cellCntr = cellCntr + 2; 931c1cad2e7SMatthew G. Knepley } 932c1cad2e7SMatthew G. Knepley } 933c1cad2e7SMatthew G. Knepley } 934c1cad2e7SMatthew G. Knepley EG_free(fobjs); 935c1cad2e7SMatthew G. Knepley } 936c1cad2e7SMatthew G. Knepley } 937c1cad2e7SMatthew G. Knepley 938c1cad2e7SMatthew G. Knepley // Generate DMPlex 9395f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, numCells, numPoints, numCorners, PETSC_TRUE, cells, cdim, coords, &dm)); 9405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(coords, cells)); 9415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(dm, " Total Number of Unique Cells = %D \n", numCells)); 9425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscInfo(dm, " Total Number of Unique Vertices = %D \n", numVertices)); 943c1cad2e7SMatthew G. Knepley 944c1cad2e7SMatthew G. Knepley // Embed EGADS model in DM 945c1cad2e7SMatthew G. Knepley { 946c1cad2e7SMatthew G. Knepley PetscContainer modelObj, contextObj; 947c1cad2e7SMatthew G. Knepley 9485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerCreate(PETSC_COMM_SELF, &modelObj)); 9495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetPointer(modelObj, model)); 9505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose((PetscObject) dm, "EGADS Model", (PetscObject) modelObj)); 9515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerDestroy(&modelObj)); 952c1cad2e7SMatthew G. Knepley 9535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerCreate(PETSC_COMM_SELF, &contextObj)); 9545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetPointer(contextObj, context)); 9555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private)); 9565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose((PetscObject) dm, "EGADS Context", (PetscObject) contextObj)); 9575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerDestroy(&contextObj)); 958c1cad2e7SMatthew G. Knepley } 959c1cad2e7SMatthew G. Knepley // Label points 960c1cad2e7SMatthew G. Knepley PetscInt nStart, nEnd; 961c1cad2e7SMatthew G. Knepley 9625f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Body ID")); 9635f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 9645f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Face ID")); 9655f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 9665f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Edge ID")); 9675f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel)); 9685f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Vertex ID")); 9695f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel)); 970c1cad2e7SMatthew G. Knepley 9715f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 2, &nStart, &nEnd)); 972c1cad2e7SMatthew G. Knepley 973c1cad2e7SMatthew G. Knepley cellCntr = 0; 974c1cad2e7SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 975c1cad2e7SMatthew G. Knepley ego body = bodies[b]; 976c1cad2e7SMatthew G. Knepley int Nv, Ne, Nf; 977c1cad2e7SMatthew G. Knepley PetscInt bodyVertexIndexStart, bodyEdgeIndexStart, bodyEdgeGlobalIndexStart, bodyFaceIndexStart; 978c1cad2e7SMatthew G. Knepley PetscHashIter BViter, BEiter, BEGiter, BFiter, EMiter; 979c1cad2e7SMatthew G. Knepley PetscBool BVfound, BEfound, BEGfound, BFfound, EMfound; 980c1cad2e7SMatthew G. Knepley 9815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyVertexMap, b, &BViter, &BVfound)); 982*28b400f6SJacob Faibussowitsch PetscCheck(BVfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyVertexMap", b); 9835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(bodyVertexMap, b, &bodyVertexIndexStart)); 984c1cad2e7SMatthew G. Knepley 9855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyEdgeMap, b, &BEiter, &BEfound)); 986*28b400f6SJacob Faibussowitsch PetscCheck(BEfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeMap", b); 9875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(bodyEdgeMap, b, &bodyEdgeIndexStart)); 988c1cad2e7SMatthew G. Knepley 9895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyFaceMap, b, &BFiter, &BFfound)); 990*28b400f6SJacob Faibussowitsch PetscCheck(BFfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyFaceMap", b); 9915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(bodyFaceMap, b, &bodyFaceIndexStart)); 992c1cad2e7SMatthew G. Knepley 9935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(bodyEdgeGlobalMap, b, &BEGiter, &BEGfound)); 994*28b400f6SJacob Faibussowitsch PetscCheck(BEGfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in bodyEdgeGlobalMap", b); 9955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(bodyEdgeGlobalMap, b, &bodyEdgeGlobalIndexStart)); 996c1cad2e7SMatthew G. Knepley 9975f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 998c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 999c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 1000c1cad2e7SMatthew G. Knepley int fID, Nl; 1001c1cad2e7SMatthew G. Knepley 10025f80ce2aSJacob Faibussowitsch fID = EG_indexBodyTopo(body, face); 1003c1cad2e7SMatthew G. Knepley 10045f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, face, LOOP, &Nl, &lobjs)); 1005c1cad2e7SMatthew G. Knepley for (int l = 0; l < Nl; ++l) { 1006c1cad2e7SMatthew G. Knepley ego loop = lobjs[l]; 1007c1cad2e7SMatthew G. Knepley int lid; 1008c1cad2e7SMatthew G. Knepley 10095f80ce2aSJacob Faibussowitsch lid = EG_indexBodyTopo(body, loop); 10102c71b3e2SJacob Faibussowitsch PetscCheckFalse(Nl > 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf); 1011c1cad2e7SMatthew G. Knepley 10125f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &eobjs, &eSenses)); 1013c1cad2e7SMatthew G. Knepley for (int e = 0; e < Ne; ++e) { 1014c1cad2e7SMatthew G. Knepley ego edge = eobjs[e]; 1015c1cad2e7SMatthew G. Knepley int eid, eOffset; 1016c1cad2e7SMatthew G. Knepley 1017c1cad2e7SMatthew G. Knepley // Skip DEGENERATE Edges 10185f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getInfo(edge, &oclass, &mtype, &topRef, &prev, &next)); 1019c1cad2e7SMatthew G. Knepley if (mtype == DEGENERATE) {continue;} 10205f80ce2aSJacob Faibussowitsch eid = EG_indexBodyTopo(body, edge); 1021c1cad2e7SMatthew G. Knepley 1022c1cad2e7SMatthew G. Knepley // get relative offset from globalEdgeID Vector 10235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &EMiter, &EMfound)); 1024*28b400f6SJacob 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); 10255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(edgeMap, bodyEdgeGlobalIndexStart + eid - 1, &eOffset)); 1026c1cad2e7SMatthew G. Knepley 10275f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, edge, NODE, &Nv, &nobjs)); 1028c1cad2e7SMatthew G. Knepley for (int v = 0; v < Nv; ++v){ 1029c1cad2e7SMatthew G. Knepley ego vertex = nobjs[v]; 1030c1cad2e7SMatthew G. Knepley int vID; 1031c1cad2e7SMatthew G. Knepley 10325f80ce2aSJacob Faibussowitsch vID = EG_indexBodyTopo(body, vertex); 10335f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, nStart + bodyVertexIndexStart + vID - 1, b)); 10345f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(vertexLabel, nStart + bodyVertexIndexStart + vID - 1, vID)); 1035c1cad2e7SMatthew G. Knepley } 1036c1cad2e7SMatthew G. Knepley EG_free(nobjs); 1037c1cad2e7SMatthew G. Knepley 10385f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, nStart + numVertices + bodyEdgeIndexStart + eOffset - 1, b)); 10395f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(edgeLabel, nStart + numVertices + bodyEdgeIndexStart + eOffset - 1, eid)); 1040c1cad2e7SMatthew G. Knepley 1041c1cad2e7SMatthew G. Knepley // Define Cell faces 1042c1cad2e7SMatthew G. Knepley for (int jj = 0; jj < 2; ++jj){ 10435f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, cellCntr, b)); 10445f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(faceLabel, cellCntr, fID)); 10455f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm, cellCntr, &cone)); 1046c1cad2e7SMatthew G. Knepley 10475f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, cone[0], b)); 10485f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(faceLabel, cone[0], fID)); 1049c1cad2e7SMatthew G. Knepley 10505f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, cone[1], b)); 10515f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(edgeLabel, cone[1], eid)); 1052c1cad2e7SMatthew G. Knepley 10535f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, cone[2], b)); 10545f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(faceLabel, cone[2], fID)); 1055c1cad2e7SMatthew G. Knepley 1056c1cad2e7SMatthew G. Knepley cellCntr = cellCntr + 1; 1057c1cad2e7SMatthew G. Knepley } 1058c1cad2e7SMatthew G. Knepley } 1059c1cad2e7SMatthew G. Knepley } 1060c1cad2e7SMatthew G. Knepley EG_free(lobjs); 1061c1cad2e7SMatthew G. Knepley 10625f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, nStart + numVertices + numEdges + bodyFaceIndexStart + fID - 1, b)); 10635f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(faceLabel, nStart + numVertices + numEdges + bodyFaceIndexStart + fID - 1, fID)); 1064c1cad2e7SMatthew G. Knepley } 1065c1cad2e7SMatthew G. Knepley EG_free(fobjs); 1066c1cad2e7SMatthew G. Knepley } 1067c1cad2e7SMatthew G. Knepley 10685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIDestroy(&edgeMap)); 10695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIDestroy(&bodyIndexMap)); 10705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIDestroy(&bodyVertexMap)); 10715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIDestroy(&bodyEdgeMap)); 10725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIDestroy(&bodyEdgeGlobalMap)); 10735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIDestroy(&bodyFaceMap)); 1074c1cad2e7SMatthew G. Knepley 1075c1cad2e7SMatthew G. Knepley *newdm = dm; 1076c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 1077c1cad2e7SMatthew G. Knepley } 1078c1cad2e7SMatthew G. Knepley 1079c1cad2e7SMatthew G. Knepley static PetscErrorCode DMPlexCreateEGADS_Tess_Internal(MPI_Comm comm, ego context, ego model, DM *newdm) 1080c1cad2e7SMatthew G. Knepley { 1081c1cad2e7SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel; 1082c1cad2e7SMatthew G. Knepley /* EGADSLite variables */ 1083c1cad2e7SMatthew G. Knepley ego geom, *bodies, *fobjs; 1084c1cad2e7SMatthew G. Knepley int b, oclass, mtype, nbodies, *senses; 1085c1cad2e7SMatthew G. Knepley int totalNumTris = 0, totalNumPoints = 0; 1086c1cad2e7SMatthew G. Knepley double boundBox[6] = {0., 0., 0., 0., 0., 0.}, tessSize; 1087c1cad2e7SMatthew G. Knepley /* PETSc variables */ 1088c1cad2e7SMatthew G. Knepley DM dm; 1089c1cad2e7SMatthew G. Knepley PetscHMapI pointIndexStartMap = NULL, triIndexStartMap = NULL, pTypeLabelMap = NULL, pIndexLabelMap = NULL; 1090c1cad2e7SMatthew G. Knepley PetscHMapI pBodyIndexLabelMap = NULL, triFaceIDLabelMap = NULL, triBodyIDLabelMap = NULL; 1091c1cad2e7SMatthew G. Knepley PetscInt dim = -1, cdim = -1, numCorners = 0, counter = 0; 1092c1cad2e7SMatthew G. Knepley PetscInt *cells = NULL; 1093c1cad2e7SMatthew G. Knepley const PetscInt *cone = NULL; 1094c1cad2e7SMatthew G. Knepley PetscReal *coords = NULL; 1095c1cad2e7SMatthew G. Knepley PetscMPIInt rank; 1096c1cad2e7SMatthew G. Knepley 1097c1cad2e7SMatthew G. Knepley PetscFunctionBeginUser; 10985f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 1099c1cad2e7SMatthew G. Knepley if (!rank) { 1100c1cad2e7SMatthew G. Knepley // --------------------------------------------------------------------------------------------------- 1101c1cad2e7SMatthew G. Knepley // Generate Petsc Plex from EGADSlite created Tessellation of geometry 1102c1cad2e7SMatthew G. Knepley // --------------------------------------------------------------------------------------------------- 1103c1cad2e7SMatthew G. Knepley 1104c1cad2e7SMatthew G. Knepley // Caluculate cell and vertex sizes 11055f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses)); 1106c1cad2e7SMatthew G. Knepley 11075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&pointIndexStartMap)); 11085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&triIndexStartMap)); 11095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&pTypeLabelMap)); 11105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&pIndexLabelMap)); 11115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&pBodyIndexLabelMap)); 11125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&triFaceIDLabelMap)); 11135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapICreate(&triBodyIDLabelMap)); 1114c1cad2e7SMatthew G. Knepley 1115c1cad2e7SMatthew G. Knepley /* Create Tessellation of Bodies */ 1116c1cad2e7SMatthew G. Knepley ego tessArray[nbodies]; 1117c1cad2e7SMatthew G. Knepley 1118c1cad2e7SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 1119c1cad2e7SMatthew G. Knepley ego body = bodies[b]; 1120c1cad2e7SMatthew G. Knepley double params[3] = {0.0, 0.0, 0.0}; // Parameters for Tessellation 1121c1cad2e7SMatthew G. Knepley int Nf, bodyNumPoints = 0, bodyNumTris = 0; 1122c1cad2e7SMatthew G. Knepley PetscHashIter PISiter, TISiter; 1123c1cad2e7SMatthew G. Knepley PetscBool PISfound, TISfound; 1124c1cad2e7SMatthew G. Knepley 1125c1cad2e7SMatthew G. Knepley /* Store Start Index for each Body's Point and Tris */ 11265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(pointIndexStartMap, b, &PISiter, &PISfound)); 11275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(triIndexStartMap, b, &TISiter, &TISfound)); 1128c1cad2e7SMatthew G. Knepley 11295f80ce2aSJacob Faibussowitsch if (!PISfound) CHKERRQ(PetscHMapISet(pointIndexStartMap, b, totalNumPoints)); 11305f80ce2aSJacob Faibussowitsch if (!TISfound) CHKERRQ(PetscHMapISet(triIndexStartMap, b, totalNumTris)); 1131c1cad2e7SMatthew G. Knepley 1132c1cad2e7SMatthew G. Knepley /* Calculate Tessellation parameters based on Bounding Box */ 1133c1cad2e7SMatthew G. Knepley /* Get Bounding Box Dimensions of the BODY */ 11345f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBoundingBox(body, boundBox)); 1135c1cad2e7SMatthew G. Knepley tessSize = boundBox[3] - boundBox[0]; 1136c1cad2e7SMatthew G. Knepley if (tessSize < boundBox[4] - boundBox[1]) tessSize = boundBox[4] - boundBox[1]; 1137c1cad2e7SMatthew G. Knepley if (tessSize < boundBox[5] - boundBox[2]) tessSize = boundBox[5] - boundBox[2]; 1138c1cad2e7SMatthew G. Knepley 1139c1cad2e7SMatthew G. Knepley // TODO :: May want to give users tessellation parameter options // 1140c1cad2e7SMatthew G. Knepley params[0] = 0.0250 * tessSize; 1141c1cad2e7SMatthew G. Knepley params[1] = 0.0075 * tessSize; 1142c1cad2e7SMatthew G. Knepley params[2] = 15.0; 1143c1cad2e7SMatthew G. Knepley 11445f80ce2aSJacob Faibussowitsch CHKERRQ(EG_makeTessBody(body, params, &tessArray[b])); 1145c1cad2e7SMatthew G. Knepley 11465f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 1147c1cad2e7SMatthew G. Knepley 1148c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 1149c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 1150c1cad2e7SMatthew G. Knepley int len, fID, ntris; 1151c1cad2e7SMatthew G. Knepley const int *ptype, *pindex, *ptris, *ptric; 1152c1cad2e7SMatthew G. Knepley const double *pxyz, *puv; 1153c1cad2e7SMatthew G. Knepley 1154c1cad2e7SMatthew G. Knepley // Get Face ID // 1155c1cad2e7SMatthew G. Knepley fID = EG_indexBodyTopo(body, face); 1156c1cad2e7SMatthew G. Knepley 1157c1cad2e7SMatthew G. Knepley // Checkout the Surface Tessellation // 11585f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric)); 1159c1cad2e7SMatthew G. Knepley 1160c1cad2e7SMatthew G. Knepley // Determine total number of triangle cells in the tessellation // 1161c1cad2e7SMatthew G. Knepley bodyNumTris += (int) ntris; 1162c1cad2e7SMatthew G. Knepley 1163c1cad2e7SMatthew G. Knepley // Check out the point index and coordinate // 1164c1cad2e7SMatthew G. Knepley for (int p = 0; p < len; ++p) { 1165c1cad2e7SMatthew G. Knepley int global; 1166c1cad2e7SMatthew G. Knepley 11675f80ce2aSJacob Faibussowitsch CHKERRQ(EG_localToGlobal(tessArray[b], fID, p+1, &global)); 1168c1cad2e7SMatthew G. Knepley 1169c1cad2e7SMatthew G. Knepley // Determine the total number of points in the tessellation // 1170c1cad2e7SMatthew G. Knepley bodyNumPoints = PetscMax(bodyNumPoints, global); 1171c1cad2e7SMatthew G. Knepley } 1172c1cad2e7SMatthew G. Knepley } 1173c1cad2e7SMatthew G. Knepley EG_free(fobjs); 1174c1cad2e7SMatthew G. Knepley 1175c1cad2e7SMatthew G. Knepley totalNumPoints += bodyNumPoints; 1176c1cad2e7SMatthew G. Knepley totalNumTris += bodyNumTris; 1177c1cad2e7SMatthew G. Knepley } 1178c1cad2e7SMatthew G. Knepley //} - Original End of (!rank) 1179c1cad2e7SMatthew G. Knepley 1180c1cad2e7SMatthew G. Knepley dim = 2; 1181c1cad2e7SMatthew G. Knepley cdim = 3; 1182c1cad2e7SMatthew G. Knepley numCorners = 3; 1183c1cad2e7SMatthew G. Knepley //PetscInt counter = 0; 1184c1cad2e7SMatthew G. Knepley 1185c1cad2e7SMatthew G. Knepley /* NEED TO DEFINE MATRICES/VECTORS TO STORE GEOM REFERENCE DATA */ 1186c1cad2e7SMatthew G. Knepley /* Fill in below and use to define DMLabels after DMPlex creation */ 11875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(totalNumPoints*cdim, &coords, totalNumTris*numCorners, &cells)); 1188c1cad2e7SMatthew G. Knepley 1189c1cad2e7SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 1190c1cad2e7SMatthew G. Knepley ego body = bodies[b]; 1191c1cad2e7SMatthew G. Knepley int Nf; 1192c1cad2e7SMatthew G. Knepley PetscInt pointIndexStart; 1193c1cad2e7SMatthew G. Knepley PetscHashIter PISiter; 1194c1cad2e7SMatthew G. Knepley PetscBool PISfound; 1195c1cad2e7SMatthew G. Knepley 11965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(pointIndexStartMap, b, &PISiter, &PISfound)); 1197*28b400f6SJacob Faibussowitsch PetscCheck(PISfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "Body %d not found in pointIndexStartMap", b); 11985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(pointIndexStartMap, b, &pointIndexStart)); 1199c1cad2e7SMatthew G. Knepley 12005f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getBodyTopos(body, NULL, FACE, &Nf, &fobjs)); 1201c1cad2e7SMatthew G. Knepley 1202c1cad2e7SMatthew G. Knepley for (int f = 0; f < Nf; ++f) { 1203c1cad2e7SMatthew G. Knepley /* Get Face Object */ 1204c1cad2e7SMatthew G. Knepley ego face = fobjs[f]; 1205c1cad2e7SMatthew G. Knepley int len, fID, ntris; 1206c1cad2e7SMatthew G. Knepley const int *ptype, *pindex, *ptris, *ptric; 1207c1cad2e7SMatthew G. Knepley const double *pxyz, *puv; 1208c1cad2e7SMatthew G. Knepley 1209c1cad2e7SMatthew G. Knepley /* Get Face ID */ 1210c1cad2e7SMatthew G. Knepley fID = EG_indexBodyTopo(body, face); 1211c1cad2e7SMatthew G. Knepley 1212c1cad2e7SMatthew G. Knepley /* Checkout the Surface Tessellation */ 12135f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTessFace(tessArray[b], fID, &len, &pxyz, &puv, &ptype, &pindex, &ntris, &ptris, &ptric)); 1214c1cad2e7SMatthew G. Knepley 1215c1cad2e7SMatthew G. Knepley /* Check out the point index and coordinate */ 1216c1cad2e7SMatthew G. Knepley for (int p = 0; p < len; ++p) { 1217c1cad2e7SMatthew G. Knepley int global; 1218c1cad2e7SMatthew G. Knepley PetscHashIter PTLiter, PILiter, PBLiter; 1219c1cad2e7SMatthew G. Knepley PetscBool PTLfound, PILfound, PBLfound; 1220c1cad2e7SMatthew G. Knepley 12215f80ce2aSJacob Faibussowitsch CHKERRQ(EG_localToGlobal(tessArray[b], fID, p+1, &global)); 1222c1cad2e7SMatthew G. Knepley 1223c1cad2e7SMatthew G. Knepley /* Set the coordinates array for DAG */ 1224c1cad2e7SMatthew G. Knepley coords[((global-1+pointIndexStart)*3) + 0] = pxyz[(p*3)+0]; 1225c1cad2e7SMatthew G. Knepley coords[((global-1+pointIndexStart)*3) + 1] = pxyz[(p*3)+1]; 1226c1cad2e7SMatthew G. Knepley coords[((global-1+pointIndexStart)*3) + 2] = pxyz[(p*3)+2]; 1227c1cad2e7SMatthew G. Knepley 1228c1cad2e7SMatthew G. Knepley /* Store Geometry Label Information for DMLabel assignment later */ 12295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(pTypeLabelMap, global-1+pointIndexStart, &PTLiter, &PTLfound)); 12305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(pIndexLabelMap, global-1+pointIndexStart, &PILiter, &PILfound)); 12315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(pBodyIndexLabelMap, global-1+pointIndexStart, &PBLiter, &PBLfound)); 1232c1cad2e7SMatthew G. Knepley 12335f80ce2aSJacob Faibussowitsch if (!PTLfound) CHKERRQ(PetscHMapISet(pTypeLabelMap, global-1+pointIndexStart, ptype[p])); 12345f80ce2aSJacob Faibussowitsch if (!PILfound) CHKERRQ(PetscHMapISet(pIndexLabelMap, global-1+pointIndexStart, pindex[p])); 12355f80ce2aSJacob Faibussowitsch if (!PBLfound) CHKERRQ(PetscHMapISet(pBodyIndexLabelMap, global-1+pointIndexStart, b)); 1236c1cad2e7SMatthew G. Knepley 12375f80ce2aSJacob Faibussowitsch if (ptype[p] < 0) CHKERRQ(PetscHMapISet(pIndexLabelMap, global-1+pointIndexStart, fID)); 1238c1cad2e7SMatthew G. Knepley } 1239c1cad2e7SMatthew G. Knepley 1240c1cad2e7SMatthew G. Knepley for (int t = 0; t < (int) ntris; ++t){ 1241c1cad2e7SMatthew G. Knepley int global, globalA, globalB; 1242c1cad2e7SMatthew G. Knepley PetscHashIter TFLiter, TBLiter; 1243c1cad2e7SMatthew G. Knepley PetscBool TFLfound, TBLfound; 1244c1cad2e7SMatthew G. Knepley 12455f80ce2aSJacob Faibussowitsch CHKERRQ(EG_localToGlobal(tessArray[b], fID, ptris[(t*3) + 0], &global)); 1246c1cad2e7SMatthew G. Knepley cells[(counter*3) +0] = global-1+pointIndexStart; 1247c1cad2e7SMatthew G. Knepley 12485f80ce2aSJacob Faibussowitsch CHKERRQ(EG_localToGlobal(tessArray[b], fID, ptris[(t*3) + 1], &globalA)); 1249c1cad2e7SMatthew G. Knepley cells[(counter*3) +1] = globalA-1+pointIndexStart; 1250c1cad2e7SMatthew G. Knepley 12515f80ce2aSJacob Faibussowitsch CHKERRQ(EG_localToGlobal(tessArray[b], fID, ptris[(t*3) + 2], &globalB)); 1252c1cad2e7SMatthew G. Knepley cells[(counter*3) +2] = globalB-1+pointIndexStart; 1253c1cad2e7SMatthew G. Knepley 12545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(triFaceIDLabelMap, counter, &TFLiter, &TFLfound)); 12555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(triBodyIDLabelMap, counter, &TBLiter, &TBLfound)); 1256c1cad2e7SMatthew G. Knepley 12575f80ce2aSJacob Faibussowitsch if (!TFLfound) CHKERRQ(PetscHMapISet(triFaceIDLabelMap, counter, fID)); 12585f80ce2aSJacob Faibussowitsch if (!TBLfound) CHKERRQ(PetscHMapISet(triBodyIDLabelMap, counter, b)); 1259c1cad2e7SMatthew G. Knepley 1260c1cad2e7SMatthew G. Knepley counter += 1; 1261c1cad2e7SMatthew G. Knepley } 1262c1cad2e7SMatthew G. Knepley } 1263c1cad2e7SMatthew G. Knepley EG_free(fobjs); 1264c1cad2e7SMatthew G. Knepley } 1265c1cad2e7SMatthew G. Knepley } 1266c1cad2e7SMatthew G. Knepley 1267c1cad2e7SMatthew G. Knepley //Build DMPlex 12685f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFromCellListPetsc(PETSC_COMM_WORLD, dim, totalNumTris, totalNumPoints, numCorners, PETSC_TRUE, cells, cdim, coords, &dm)); 12695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(coords, cells)); 1270c1cad2e7SMatthew G. Knepley 1271c1cad2e7SMatthew G. Knepley // Embed EGADS model in DM 1272c1cad2e7SMatthew G. Knepley { 1273c1cad2e7SMatthew G. Knepley PetscContainer modelObj, contextObj; 1274c1cad2e7SMatthew G. Knepley 12755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerCreate(PETSC_COMM_SELF, &modelObj)); 12765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetPointer(modelObj, model)); 12775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose((PetscObject) dm, "EGADS Model", (PetscObject) modelObj)); 12785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerDestroy(&modelObj)); 1279c1cad2e7SMatthew G. Knepley 12805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerCreate(PETSC_COMM_SELF, &contextObj)); 12815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetPointer(contextObj, context)); 12825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerSetUserDestroy(contextObj, DMPlexEGADSDestroy_Private)); 12835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectCompose((PetscObject) dm, "EGADS Context", (PetscObject) contextObj)); 12845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerDestroy(&contextObj)); 1285c1cad2e7SMatthew G. Knepley } 1286c1cad2e7SMatthew G. Knepley 1287c1cad2e7SMatthew G. Knepley // Label Points 12885f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Body ID")); 12895f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 12905f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Face ID")); 12915f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 12925f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Edge ID")); 12935f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel)); 12945f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(dm, "EGADS Vertex ID")); 12955f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel)); 1296c1cad2e7SMatthew G. Knepley 1297c1cad2e7SMatthew G. Knepley /* Get Number of DAG Nodes at each level */ 1298c1cad2e7SMatthew G. Knepley int fStart, fEnd, eStart, eEnd, nStart, nEnd; 1299c1cad2e7SMatthew G. Knepley 13005f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 0, &fStart, &fEnd)); 13015f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 1, &eStart, &eEnd)); 13025f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetHeightStratum(dm, 2, &nStart, &nEnd)); 1303c1cad2e7SMatthew G. Knepley 1304c1cad2e7SMatthew G. Knepley /* Set DMLabels for NODES */ 1305c1cad2e7SMatthew G. Knepley for (int n = nStart; n < nEnd; ++n) { 1306c1cad2e7SMatthew G. Knepley int pTypeVal, pIndexVal, pBodyVal; 1307c1cad2e7SMatthew G. Knepley PetscHashIter PTLiter, PILiter, PBLiter; 1308c1cad2e7SMatthew G. Knepley PetscBool PTLfound, PILfound, PBLfound; 1309c1cad2e7SMatthew G. Knepley 1310c1cad2e7SMatthew G. Knepley //Converted to Hash Tables 13115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(pTypeLabelMap, n - nStart, &PTLiter, &PTLfound)); 1312*28b400f6SJacob Faibussowitsch PetscCheck(PTLfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pTypeLabelMap", n); 13135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(pTypeLabelMap, n - nStart, &pTypeVal)); 1314c1cad2e7SMatthew G. Knepley 13155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(pIndexLabelMap, n - nStart, &PILiter, &PILfound)); 1316*28b400f6SJacob Faibussowitsch PetscCheck(PILfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pIndexLabelMap", n); 13175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(pIndexLabelMap, n - nStart, &pIndexVal)); 1318c1cad2e7SMatthew G. Knepley 13195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(pBodyIndexLabelMap, n - nStart, &PBLiter, &PBLfound)); 1320*28b400f6SJacob Faibussowitsch PetscCheck(PBLfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in pBodyLabelMap", n); 13215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(pBodyIndexLabelMap, n - nStart, &pBodyVal)); 1322c1cad2e7SMatthew G. Knepley 13235f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, n, pBodyVal)); 13245f80ce2aSJacob Faibussowitsch if (pTypeVal == 0) CHKERRQ(DMLabelSetValue(vertexLabel, n, pIndexVal)); 13255f80ce2aSJacob Faibussowitsch if (pTypeVal > 0) CHKERRQ(DMLabelSetValue(edgeLabel, n, pIndexVal)); 13265f80ce2aSJacob Faibussowitsch if (pTypeVal < 0) CHKERRQ(DMLabelSetValue(faceLabel, n, pIndexVal)); 1327c1cad2e7SMatthew G. Knepley } 1328c1cad2e7SMatthew G. Knepley 1329c1cad2e7SMatthew G. Knepley /* Set DMLabels for Edges - Based on the DMLabels of the EDGE's NODES */ 1330c1cad2e7SMatthew G. Knepley for (int e = eStart; e < eEnd; ++e) { 1331c1cad2e7SMatthew G. Knepley int bodyID_0, vertexID_0, vertexID_1, edgeID_0, edgeID_1, faceID_0, faceID_1; 1332c1cad2e7SMatthew G. Knepley 13335f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm, e, &cone)); 13345f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(bodyLabel, cone[0], &bodyID_0)); // Do I need to check the other end? 13355f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(vertexLabel, cone[0], &vertexID_0)); 13365f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(vertexLabel, cone[1], &vertexID_1)); 13375f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(edgeLabel, cone[0], &edgeID_0)); 13385f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(edgeLabel, cone[1], &edgeID_1)); 13395f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(faceLabel, cone[0], &faceID_0)); 13405f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(faceLabel, cone[1], &faceID_1)); 1341c1cad2e7SMatthew G. Knepley 13425f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, e, bodyID_0)); 1343c1cad2e7SMatthew G. Knepley 13445f80ce2aSJacob Faibussowitsch if (edgeID_0 == edgeID_1) CHKERRQ(DMLabelSetValue(edgeLabel, e, edgeID_0)); 13455f80ce2aSJacob Faibussowitsch else if (vertexID_0 > 0 && edgeID_1 > 0) CHKERRQ(DMLabelSetValue(edgeLabel, e, edgeID_1)); 13465f80ce2aSJacob Faibussowitsch else if (vertexID_1 > 0 && edgeID_0 > 0) CHKERRQ(DMLabelSetValue(edgeLabel, e, edgeID_0)); 1347c1cad2e7SMatthew G. Knepley else { /* Do Nothing */ } 1348c1cad2e7SMatthew G. Knepley } 1349c1cad2e7SMatthew G. Knepley 1350c1cad2e7SMatthew G. Knepley /* Set DMLabels for Cells */ 1351c1cad2e7SMatthew G. Knepley for (int f = fStart; f < fEnd; ++f){ 1352c1cad2e7SMatthew G. Knepley int edgeID_0; 1353c1cad2e7SMatthew G. Knepley PetscInt triBodyVal, triFaceVal; 1354c1cad2e7SMatthew G. Knepley PetscHashIter TFLiter, TBLiter; 1355c1cad2e7SMatthew G. Knepley PetscBool TFLfound, TBLfound; 1356c1cad2e7SMatthew G. Knepley 1357c1cad2e7SMatthew G. Knepley // Convert to Hash Table 13585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(triFaceIDLabelMap, f - fStart, &TFLiter, &TFLfound)); 1359*28b400f6SJacob Faibussowitsch PetscCheck(TFLfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in triFaceIDLabelMap", f); 13605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(triFaceIDLabelMap, f - fStart, &triFaceVal)); 1361c1cad2e7SMatthew G. Knepley 13625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIFind(triBodyIDLabelMap, f - fStart, &TBLiter, &TBLfound)); 1363*28b400f6SJacob Faibussowitsch PetscCheck(TBLfound,PETSC_COMM_SELF, PETSC_ERR_SUP, "DAG Point %d not found in triBodyIDLabelMap", f); 13645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscHMapIGet(triBodyIDLabelMap, f - fStart, &triBodyVal)); 1365c1cad2e7SMatthew G. Knepley 13665f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, f, triBodyVal)); 13675f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(faceLabel, f, triFaceVal)); 1368c1cad2e7SMatthew G. Knepley 1369c1cad2e7SMatthew G. Knepley /* Finish Labeling previously unlabeled DMPlex Edges - Assumes Triangular Cell (3 Edges Max) */ 13705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetCone(dm, f, &cone)); 1371c1cad2e7SMatthew G. Knepley 1372c1cad2e7SMatthew G. Knepley for (int jj = 0; jj < 3; ++jj) { 13735f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(edgeLabel, cone[jj], &edgeID_0)); 1374c1cad2e7SMatthew G. Knepley 1375c1cad2e7SMatthew G. Knepley if (edgeID_0 < 0) { 13765f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(bodyLabel, cone[jj], triBodyVal)); 13775f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(faceLabel, cone[jj], triFaceVal)); 1378c1cad2e7SMatthew G. Knepley } 1379c1cad2e7SMatthew G. Knepley } 1380c1cad2e7SMatthew G. Knepley } 1381c1cad2e7SMatthew G. Knepley 1382c1cad2e7SMatthew G. Knepley *newdm = dm; 1383c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 1384c1cad2e7SMatthew G. Knepley } 13857bee2925SMatthew Knepley #endif 13867bee2925SMatthew Knepley 1387c1cad2e7SMatthew G. Knepley /*@ 1388c1cad2e7SMatthew 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. 1389c1cad2e7SMatthew G. Knepley 1390c1cad2e7SMatthew G. Knepley Collective on dm 1391c1cad2e7SMatthew G. Knepley 1392c1cad2e7SMatthew G. Knepley Input Parameter: 1393c1cad2e7SMatthew G. Knepley . dm - The uninflated DM object representing the mesh 1394c1cad2e7SMatthew G. Knepley 1395c1cad2e7SMatthew G. Knepley Output Parameter: 1396c1cad2e7SMatthew G. Knepley . dm - The inflated DM object representing the mesh 1397c1cad2e7SMatthew G. Knepley 1398c1cad2e7SMatthew G. Knepley Level: intermediate 1399c1cad2e7SMatthew G. Knepley 1400c1cad2e7SMatthew G. Knepley .seealso: DMPLEX, DMCreate(), DMPlexCreateEGADS() 1401c1cad2e7SMatthew G. Knepley @*/ 1402c1cad2e7SMatthew G. Knepley PetscErrorCode DMPlexInflateToGeomModel(DM dm) 1403c1cad2e7SMatthew G. Knepley { 1404c1cad2e7SMatthew G. Knepley #if defined(PETSC_HAVE_EGADS) 1405c1cad2e7SMatthew G. Knepley /* EGADS Variables */ 1406c1cad2e7SMatthew G. Knepley ego model, geom, body, face, edge; 1407c1cad2e7SMatthew G. Knepley ego *bodies; 1408c1cad2e7SMatthew G. Knepley int Nb, oclass, mtype, *senses; 1409c1cad2e7SMatthew G. Knepley double result[3]; 1410c1cad2e7SMatthew G. Knepley /* PETSc Variables */ 1411c1cad2e7SMatthew G. Knepley DM cdm; 1412c1cad2e7SMatthew G. Knepley PetscContainer modelObj; 1413c1cad2e7SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel, vertexLabel; 1414c1cad2e7SMatthew G. Knepley Vec coordinates; 1415c1cad2e7SMatthew G. Knepley PetscScalar *coords; 1416c1cad2e7SMatthew G. Knepley PetscInt bodyID, faceID, edgeID, vertexID; 1417c1cad2e7SMatthew G. Knepley PetscInt cdim, d, vStart, vEnd, v; 1418c1cad2e7SMatthew G. Knepley #endif 1419c1cad2e7SMatthew G. Knepley 1420c1cad2e7SMatthew G. Knepley PetscFunctionBegin; 1421c1cad2e7SMatthew G. Knepley #if defined(PETSC_HAVE_EGADS) 14225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectQuery((PetscObject) dm, "EGADS Model", (PetscObject *) &modelObj)); 1423c1cad2e7SMatthew G. Knepley if (!modelObj) PetscFunctionReturn(0); 14245f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDim(dm, &cdim)); 14255f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(dm, &cdm)); 14265f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(dm, &coordinates)); 14275f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Body ID", &bodyLabel)); 14285f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Face ID", &faceLabel)); 14295f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Edge ID", &edgeLabel)); 14305f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(dm, "EGADS Vertex ID", &vertexLabel)); 1431c1cad2e7SMatthew G. Knepley 14325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscContainerGetPointer(modelObj, (void **) &model)); 14335f80ce2aSJacob Faibussowitsch CHKERRQ(EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses)); 1434c1cad2e7SMatthew G. Knepley 14355f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 14365f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(coordinates, &coords)); 1437c1cad2e7SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 1438c1cad2e7SMatthew G. Knepley PetscScalar *vcoords; 1439c1cad2e7SMatthew G. Knepley 14405f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(bodyLabel, v, &bodyID)); 14415f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(faceLabel, v, &faceID)); 14425f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(edgeLabel, v, &edgeID)); 14435f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelGetValue(vertexLabel, v, &vertexID)); 1444c1cad2e7SMatthew G. Knepley 14452c71b3e2SJacob Faibussowitsch PetscCheckFalse(bodyID >= Nb,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Body %D is not in [0, %d)", bodyID, Nb); 1446c1cad2e7SMatthew G. Knepley body = bodies[bodyID]; 1447c1cad2e7SMatthew G. Knepley 14485f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexPointLocalRef(cdm, v, coords, (void *) &vcoords)); 1449c1cad2e7SMatthew G. Knepley if (edgeID > 0) { 1450c1cad2e7SMatthew G. Knepley /* Snap to EDGE at nearest location */ 1451c1cad2e7SMatthew G. Knepley double params[1]; 14525f80ce2aSJacob Faibussowitsch CHKERRQ(EG_objectBodyTopo(body, EDGE, edgeID, &edge)); 14535f80ce2aSJacob Faibussowitsch CHKERRQ(EG_invEvaluate(edge, vcoords, params, result)); // Get (x,y,z) of nearest point on EDGE 1454c1cad2e7SMatthew G. Knepley for (d = 0; d < cdim; ++d) vcoords[d] = result[d]; 1455c1cad2e7SMatthew G. Knepley } else if (faceID > 0) { 1456c1cad2e7SMatthew G. Knepley /* Snap to FACE at nearest location */ 1457c1cad2e7SMatthew G. Knepley double params[2]; 14585f80ce2aSJacob Faibussowitsch CHKERRQ(EG_objectBodyTopo(body, FACE, faceID, &face)); 14595f80ce2aSJacob Faibussowitsch CHKERRQ(EG_invEvaluate(face, vcoords, params, result)); // Get (x,y,z) of nearest point on FACE 1460c1cad2e7SMatthew G. Knepley for (d = 0; d < cdim; ++d) vcoords[d] = result[d]; 1461c1cad2e7SMatthew G. Knepley } 1462c1cad2e7SMatthew G. Knepley } 14635f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(coordinates, &coords)); 1464c1cad2e7SMatthew G. Knepley /* Clear out global coordinates */ 14655f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&dm->coordinates)); 1466c1cad2e7SMatthew G. Knepley #endif 1467c1cad2e7SMatthew G. Knepley PetscFunctionReturn(0); 1468c1cad2e7SMatthew G. Knepley } 1469c1cad2e7SMatthew G. Knepley 14707bee2925SMatthew Knepley /*@C 1471c1cad2e7SMatthew G. Knepley DMPlexCreateEGADSFromFile - Create a DMPlex mesh from an EGADS, IGES, or STEP file. 14727bee2925SMatthew Knepley 14737bee2925SMatthew Knepley Collective 14747bee2925SMatthew Knepley 14757bee2925SMatthew Knepley Input Parameters: 14767bee2925SMatthew Knepley + comm - The MPI communicator 1477c1cad2e7SMatthew G. Knepley - filename - The name of the EGADS, IGES, or STEP file 14787bee2925SMatthew Knepley 14797bee2925SMatthew Knepley Output Parameter: 14807bee2925SMatthew Knepley . dm - The DM object representing the mesh 14817bee2925SMatthew Knepley 14827bee2925SMatthew Knepley Level: beginner 14837bee2925SMatthew Knepley 1484c1cad2e7SMatthew G. Knepley .seealso: DMPLEX, DMCreate(), DMPlexCreateEGADS(), DMPlexCreateEGADSLiteFromFile() 14857bee2925SMatthew Knepley @*/ 14867bee2925SMatthew Knepley PetscErrorCode DMPlexCreateEGADSFromFile(MPI_Comm comm, const char filename[], DM *dm) 14877bee2925SMatthew Knepley { 14887bee2925SMatthew Knepley PetscMPIInt rank; 14897bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS) 14907bee2925SMatthew Knepley ego context= NULL, model = NULL; 14917bee2925SMatthew Knepley #endif 1492c1cad2e7SMatthew G. Knepley PetscBool printModel = PETSC_FALSE, tessModel = PETSC_FALSE, newModel = PETSC_FALSE; 14937bee2925SMatthew Knepley 14947bee2925SMatthew Knepley PetscFunctionBegin; 14957bee2925SMatthew Knepley PetscValidCharPointer(filename, 2); 14965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_print_model", &printModel, NULL)); 14975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_tess_model", &tessModel, NULL)); 14985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL, NULL, "-dm_plex_egads_new_model", &newModel, NULL)); 14995f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 15007bee2925SMatthew Knepley #if defined(PETSC_HAVE_EGADS) 1501dd400576SPatrick Sanan if (rank == 0) { 15027bee2925SMatthew Knepley 15035f80ce2aSJacob Faibussowitsch CHKERRQ(EG_open(&context)); 15045f80ce2aSJacob Faibussowitsch CHKERRQ(EG_loadModel(context, 0, filename, &model)); 15055f80ce2aSJacob Faibussowitsch if (printModel) CHKERRQ(DMPlexEGADSPrintModel_Internal(model)); 15067bee2925SMatthew Knepley 15077bee2925SMatthew Knepley } 15085f80ce2aSJacob Faibussowitsch if (tessModel) CHKERRQ(DMPlexCreateEGADS_Tess_Internal(comm, context, model, dm)); 15095f80ce2aSJacob Faibussowitsch else if (newModel) CHKERRQ(DMPlexCreateEGADS(comm, context, model, dm)); 15105f80ce2aSJacob Faibussowitsch else CHKERRQ(DMPlexCreateEGADS_Internal(comm, context, model, dm)); 1511c03d1177SJose E. Roman PetscFunctionReturn(0); 15127bee2925SMatthew Knepley #else 1513c1cad2e7SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "This method requires EGADS support. Reconfigure using --download-egads"); 15147bee2925SMatthew Knepley #endif 15157bee2925SMatthew Knepley } 1516