1c4762a1bSJed Brown static const char help[] = "Test of EGADSLite CAD functionality"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscdmplex.h> 46f4f5c14SMatthew G. Knepley #include <petsc/private/hashmapi.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown #include <egads.h> 7c4762a1bSJed Brown #include <petsc.h> 8c4762a1bSJed Brown 9c4762a1bSJed Brown typedef struct { 10c4762a1bSJed Brown char filename[PETSC_MAX_PATH_LEN]; 11c4762a1bSJed Brown } AppCtx; 12c4762a1bSJed Brown 13c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 14c4762a1bSJed Brown { 15c4762a1bSJed Brown PetscErrorCode ierr; 16c4762a1bSJed Brown 17c4762a1bSJed Brown PetscFunctionBeginUser; 18c4762a1bSJed Brown options->filename[0] = '\0'; 19c4762a1bSJed Brown 20c4762a1bSJed Brown ierr = PetscOptionsBegin(comm, "", "EGADSPlex Problem Options", "EGADSLite");CHKERRQ(ierr); 21*589a23caSBarry Smith ierr = PetscOptionsString("-filename", "The EGADSLite file", "ex9.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr); 22c4762a1bSJed Brown ierr = PetscOptionsEnd(); 23c4762a1bSJed Brown PetscFunctionReturn(0); 24c4762a1bSJed Brown } 25c4762a1bSJed Brown 266f4f5c14SMatthew G. Knepley static PetscErrorCode PrintEGADS(ego model) 27c4762a1bSJed Brown { 286f4f5c14SMatthew G. Knepley ego geom, *bodies, *objs, *nobjs, *mobjs, *lobjs; 296f4f5c14SMatthew G. Knepley int oclass, mtype, *senses; 306f4f5c14SMatthew G. Knepley int Nb, b; 31c4762a1bSJed Brown PetscErrorCode ierr; 32c4762a1bSJed Brown 336f4f5c14SMatthew G. Knepley PetscFunctionBeginUser; 34c4762a1bSJed Brown /* test bodyTopo functions */ 356f4f5c14SMatthew G. Knepley ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &Nb, &bodies, &senses);CHKERRQ(ierr); 366f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " Number of BODIES (nbodies): %d \n", Nb);CHKERRQ(ierr); 37c4762a1bSJed Brown 386f4f5c14SMatthew G. Knepley for (b = 0; b < Nb; ++b) { 39c4762a1bSJed Brown ego body = bodies[b]; 40c4762a1bSJed Brown int id, Nsh, Nf, Nl, l, Ne, e, Nv, v; 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* Output Basic Model Topology */ 43c4762a1bSJed Brown ierr = EG_getBodyTopos(body, NULL, SHELL, &Nsh, &objs);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " Number of SHELLS: %d \n", Nsh);CHKERRQ(ierr); 456f4f5c14SMatthew G. Knepley EG_free(objs); 46c4762a1bSJed Brown 47c4762a1bSJed Brown ierr = EG_getBodyTopos(body, NULL, FACE, &Nf, &objs);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " Number of FACES: %d \n", Nf);CHKERRQ(ierr); 496f4f5c14SMatthew G. Knepley EG_free(objs); 50c4762a1bSJed Brown 51c4762a1bSJed Brown ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " Number of LOOPS: %d \n", Nl);CHKERRQ(ierr); 53c4762a1bSJed Brown 54c4762a1bSJed Brown ierr = EG_getBodyTopos(body, NULL, EDGE, &Ne, &objs);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " Number of EDGES: %d \n", Ne);CHKERRQ(ierr); 566f4f5c14SMatthew G. Knepley EG_free(objs); 57c4762a1bSJed Brown 58c4762a1bSJed Brown ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &objs);CHKERRQ(ierr); 59c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " Number of NODES: %d \n", Nv);CHKERRQ(ierr); 606f4f5c14SMatthew G. Knepley EG_free(objs); 61c4762a1bSJed Brown 62c4762a1bSJed Brown for (l = 0; l < Nl; ++l) { 63c4762a1bSJed Brown ego loop = lobjs[l]; 64c4762a1bSJed Brown 65c4762a1bSJed Brown id = EG_indexBodyTopo(body, loop); 66c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " LOOP ID: %d\n", id);CHKERRQ(ierr); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* Get EDGE info which associated with the current LOOP */ 69c4762a1bSJed Brown ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr); 70c4762a1bSJed Brown 71c4762a1bSJed Brown for (e = 0; e < Ne; ++e) { 72c4762a1bSJed Brown ego edge = objs[e]; 73c4762a1bSJed Brown 74c4762a1bSJed Brown id = EG_indexBodyTopo(body, edge);CHKERRQ(ierr); 75c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " EDGE ID: %d\n", id);CHKERRQ(ierr); 76c4762a1bSJed Brown 77c4762a1bSJed Brown double range[4] = {0., 0., 0., 0.}; 78c4762a1bSJed Brown double point[3] = {0., 0., 0.}; 79c4762a1bSJed Brown int peri; 80c4762a1bSJed Brown 81c4762a1bSJed Brown ierr = EG_getRange(objs[e], range, &peri); 82c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " Range = %lf, %lf, %lf, %lf \n", range[0], range[1], range[2], range[3]); 83c4762a1bSJed Brown 84c4762a1bSJed Brown /* Get NODE info which associated with the current EDGE */ 85c4762a1bSJed Brown ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses);CHKERRQ(ierr); 86c4762a1bSJed Brown 87c4762a1bSJed Brown for (v = 0; v < Nv; ++v) { 88c4762a1bSJed Brown ego vertex = nobjs[v]; 89c4762a1bSJed Brown double limits[4]; 90c4762a1bSJed Brown int dummy; 91c4762a1bSJed Brown 92c4762a1bSJed Brown ierr = EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses);CHKERRQ(ierr); 93c4762a1bSJed Brown id = EG_indexBodyTopo(body, vertex); 94c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " NODE ID: %d \n", id);CHKERRQ(ierr); 95c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " (x, y, z) = (%lf, %lf, %lf) \n", limits[0], limits[1], limits[2]); 96c4762a1bSJed Brown 97c4762a1bSJed Brown point[0] = point[0] + limits[0]; 98c4762a1bSJed Brown point[1] = point[1] + limits[1]; 99c4762a1bSJed Brown point[2] = point[2] + limits[2]; 100c4762a1bSJed Brown } 101c4762a1bSJed Brown } 102c4762a1bSJed Brown } 1036f4f5c14SMatthew G. Knepley EG_free(lobjs); 104c4762a1bSJed Brown } 1056f4f5c14SMatthew G. Knepley PetscFunctionReturn(0); 1066f4f5c14SMatthew G. Knepley } 1076f4f5c14SMatthew G. Knepley 1086f4f5c14SMatthew G. Knepley int main(int argc, char *argv[]) 1096f4f5c14SMatthew G. Knepley { 1106f4f5c14SMatthew G. Knepley DMLabel bodyLabel, faceLabel, edgeLabel; 1116f4f5c14SMatthew G. Knepley PetscInt cStart, cEnd, c; 1126f4f5c14SMatthew G. Knepley /* EGADSLite variables */ 1136f4f5c14SMatthew G. Knepley ego context, model, geom, *bodies, *objs, *nobjs, *mobjs, *lobjs; 1146f4f5c14SMatthew G. Knepley int oclass, mtype, nbodies, *senses; 1156f4f5c14SMatthew G. Knepley int b; 1166f4f5c14SMatthew G. Knepley /* PETSc variables */ 1176f4f5c14SMatthew G. Knepley DM dm; 1186f4f5c14SMatthew G. Knepley PetscHMapI edgeMap = NULL; 1196f4f5c14SMatthew G. 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; 1206f4f5c14SMatthew G. Knepley PetscInt *cells = NULL, *cone = NULL; 1216f4f5c14SMatthew G. Knepley PetscReal *coords = NULL; 1226f4f5c14SMatthew G. Knepley MPI_Comm comm; 1236f4f5c14SMatthew G. Knepley PetscMPIInt rank; 1246f4f5c14SMatthew G. Knepley AppCtx ctx; 1256f4f5c14SMatthew G. Knepley PetscErrorCode ierr; 1266f4f5c14SMatthew G. Knepley 1276f4f5c14SMatthew G. Knepley ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr; 1286f4f5c14SMatthew G. Knepley comm = PETSC_COMM_WORLD; 1296f4f5c14SMatthew G. Knepley ierr = ProcessOptions(comm, &ctx);CHKERRQ(ierr); 1306f4f5c14SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1316f4f5c14SMatthew G. Knepley if (!rank) { 1326f4f5c14SMatthew G. Knepley /* Open EGADs file and load EGADs model data */ 1336f4f5c14SMatthew G. Knepley ierr = EG_open(&context);CHKERRQ(ierr); 1346f4f5c14SMatthew G. Knepley ierr = EG_loadModel(context, 0, ctx.filename, &model);CHKERRQ(ierr); 1356f4f5c14SMatthew G. Knepley 1366f4f5c14SMatthew G. Knepley ierr = PrintEGADS(model);CHKERRQ(ierr); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* --------------------------------------------------------------------------------------------------- 139c4762a1bSJed Brown Generate Petsc Plex 140c4762a1bSJed Brown Get all Nodes in model, record coordinates in a correctly formatted array 1416f4f5c14SMatthew G. Knepley Cycle through bodies, cycle through loops, recorde NODE IDs in a correctly formatted array 1426f4f5c14SMatthew G. Knepley We need to uniformly refine the initial geometry to guarantee a valid mesh 1436f4f5c14SMatthew G. Knepley */ 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Calculate cell and vertex sizes */ 146c4762a1bSJed Brown ierr = EG_getTopology(model, &geom, &oclass, &mtype, NULL, &nbodies, &bodies, &senses);CHKERRQ(ierr); 1476f4f5c14SMatthew G. Knepley ierr = PetscHMapICreate(&edgeMap);CHKERRQ(ierr); 1486f4f5c14SMatthew G. Knepley numEdges = 0; 149c4762a1bSJed Brown for (b = 0; b < nbodies; ++b) { 150c4762a1bSJed Brown ego body = bodies[b]; 151c4762a1bSJed Brown int id, Nl, l, Nv, v; 152c4762a1bSJed Brown 153c4762a1bSJed Brown ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr); 154c4762a1bSJed Brown for (l = 0; l < Nl; ++l) { 155c4762a1bSJed Brown ego loop = lobjs[l]; 1566f4f5c14SMatthew G. Knepley int Ner = 0, Ne, e, Nc; 157c4762a1bSJed Brown 1586f4f5c14SMatthew G. Knepley ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr); 1596f4f5c14SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 1606f4f5c14SMatthew G. Knepley ego edge = objs[e]; 1616f4f5c14SMatthew G. Knepley int Nv, id; 1626f4f5c14SMatthew G. Knepley PetscHashIter iter; 1636f4f5c14SMatthew G. Knepley PetscBool found; 1646f4f5c14SMatthew G. Knepley 1656f4f5c14SMatthew G. Knepley ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses); 1666f4f5c14SMatthew G. Knepley if (mtype == DEGENERATE) continue; 1676f4f5c14SMatthew G. Knepley id = EG_indexBodyTopo(body, edge);CHKERRQ(ierr); 1686f4f5c14SMatthew G. Knepley ierr = PetscHMapIFind(edgeMap, id-1, &iter, &found);CHKERRQ(ierr); 1696f4f5c14SMatthew G. Knepley if (!found) {ierr = PetscHMapISet(edgeMap, id-1, numEdges++);CHKERRQ(ierr);} 1706f4f5c14SMatthew G. Knepley ++Ner; 171c4762a1bSJed Brown } 1726f4f5c14SMatthew G. Knepley if (Ner == 2) {Nc = 2;} 1736f4f5c14SMatthew G. Knepley else if (Ner == 3) {Nc = 4;} 1746f4f5c14SMatthew G. Knepley else if (Ner == 4) {Nc = 8; ++numQuads;} 1756f4f5c14SMatthew G. Knepley else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot support loop with %d edges", Ner); 1766f4f5c14SMatthew G. Knepley numCells += Nc; 1776f4f5c14SMatthew G. Knepley newCells += Nc-1; 1786f4f5c14SMatthew G. Knepley maxCorners = PetscMax(Ner*2+1, maxCorners); 1796f4f5c14SMatthew G. Knepley } 1806f4f5c14SMatthew G. Knepley ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs);CHKERRQ(ierr); 181c4762a1bSJed Brown for (v = 0; v < Nv; ++v) { 182c4762a1bSJed Brown ego vertex = nobjs[v]; 183c4762a1bSJed Brown 184c4762a1bSJed Brown id = EG_indexBodyTopo(body, vertex); 185c4762a1bSJed Brown /* TODO: Instead of assuming contiguous ids, we could use a hash table */ 186c4762a1bSJed Brown numVertices = PetscMax(id, numVertices); 187c4762a1bSJed Brown } 1886f4f5c14SMatthew G. Knepley EG_free(lobjs); 1896f4f5c14SMatthew G. Knepley EG_free(nobjs); 190c4762a1bSJed Brown } 1916f4f5c14SMatthew G. Knepley ierr = PetscHMapIGetSize(edgeMap, &numEdges);CHKERRQ(ierr); 1926f4f5c14SMatthew G. Knepley newVertices = numEdges + numQuads; 1936f4f5c14SMatthew G. Knepley numVertices += newVertices; 194c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, "\nPLEX Input Array Checkouts\n");CHKERRQ(ierr); 1956f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " Total Number of Unique Cells = %D (%D)\n", numCells, newCells);CHKERRQ(ierr); 1966f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " Total Number of Unique Vertices = %D (%D) \n", numVertices, newVertices);CHKERRQ(ierr); 197c4762a1bSJed Brown 198c4762a1bSJed Brown dim = 2; /* Assume 3D Models :: Need to update to handle 2D Models in the future */ 199c4762a1bSJed Brown cdim = 3; /* Assume 3D Models :: Need to update to handle 2D Models in the future */ 2006f4f5c14SMatthew G. Knepley numCorners = 3; /* Split cells into triangles */ 2016f4f5c14SMatthew G. Knepley ierr = PetscMalloc3(numVertices*cdim, &coords, numCells*numCorners, &cells, maxCorners, &cone);CHKERRQ(ierr); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* Get vertex coordinates */ 204c4762a1bSJed Brown for (b = 0; b < nbodies; ++b) { 205c4762a1bSJed Brown ego body = bodies[b]; 206c4762a1bSJed Brown int id, Nv, v; 207c4762a1bSJed Brown 208c4762a1bSJed Brown ierr = EG_getBodyTopos(body, NULL, NODE, &Nv, &nobjs);CHKERRQ(ierr); 209c4762a1bSJed Brown for (v = 0; v < Nv; ++v) { 210c4762a1bSJed Brown ego vertex = nobjs[v]; 211c4762a1bSJed Brown double limits[4]; 212c4762a1bSJed Brown int dummy; 213c4762a1bSJed Brown 214c4762a1bSJed Brown ierr = EG_getTopology(vertex, &geom, &oclass, &mtype, limits, &dummy, &mobjs, &senses);CHKERRQ(ierr); 215c4762a1bSJed Brown id = EG_indexBodyTopo(body, vertex);CHKERRQ(ierr); 216c4762a1bSJed Brown coords[(id-1)*cdim+0] = limits[0]; 217c4762a1bSJed Brown coords[(id-1)*cdim+1] = limits[1]; 218c4762a1bSJed Brown coords[(id-1)*cdim+2] = limits[2]; 2196f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " Node ID = %d (%d)\n", id, id-1); 220c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " (x,y,z) = (%lf, %lf, %lf) \n \n", coords[(id-1)*cdim+0], coords[(id-1)*cdim+1],coords[(id-1)*cdim+2]); 221c4762a1bSJed Brown } 2226f4f5c14SMatthew G. Knepley EG_free(nobjs); 223c4762a1bSJed Brown } 2246f4f5c14SMatthew G. Knepley ierr = PetscHMapIClear(edgeMap);CHKERRQ(ierr); 2256f4f5c14SMatthew G. Knepley fOff = numVertices - newVertices + numEdges; 2266f4f5c14SMatthew G. Knepley numEdges = 0; 2276f4f5c14SMatthew G. Knepley numQuads = 0; 2286f4f5c14SMatthew G. Knepley for (b = 0; b < nbodies; ++b) { 2296f4f5c14SMatthew G. Knepley ego body = bodies[b]; 2306f4f5c14SMatthew G. Knepley int Nl, l; 2316f4f5c14SMatthew G. Knepley 2326f4f5c14SMatthew G. Knepley ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr); 2336f4f5c14SMatthew G. Knepley for (l = 0; l < Nl; ++l) { 2346f4f5c14SMatthew G. Knepley ego loop = lobjs[l]; 2356f4f5c14SMatthew G. Knepley int lid, Ner = 0, Ne, e; 2366f4f5c14SMatthew G. Knepley 2376f4f5c14SMatthew G. Knepley lid = EG_indexBodyTopo(body, loop);CHKERRQ(ierr); 2386f4f5c14SMatthew G. Knepley ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr); 2396f4f5c14SMatthew G. Knepley for (e = 0; e < Ne; ++e) { 2406f4f5c14SMatthew G. Knepley ego edge = objs[e]; 2416f4f5c14SMatthew G. Knepley int eid, Nv; 2426f4f5c14SMatthew G. Knepley PetscHashIter iter; 2436f4f5c14SMatthew G. Knepley PetscBool found; 2446f4f5c14SMatthew G. Knepley 2456f4f5c14SMatthew G. Knepley ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses); 2466f4f5c14SMatthew G. Knepley if (mtype == DEGENERATE) continue; 2476f4f5c14SMatthew G. Knepley ++Ner; 2486f4f5c14SMatthew G. Knepley eid = EG_indexBodyTopo(body, edge);CHKERRQ(ierr); 2496f4f5c14SMatthew G. Knepley ierr = PetscHMapIFind(edgeMap, eid-1, &iter, &found);CHKERRQ(ierr); 2506f4f5c14SMatthew G. Knepley if (!found) { 2516f4f5c14SMatthew G. Knepley PetscInt v = numVertices - newVertices + numEdges; 2526f4f5c14SMatthew G. Knepley double range[4], params[3] = {0., 0., 0.}, result[18]; 2536f4f5c14SMatthew G. Knepley int periodic[2]; 2546f4f5c14SMatthew G. Knepley 2556f4f5c14SMatthew G. Knepley ierr = PetscHMapISet(edgeMap, eid-1, numEdges++);CHKERRQ(ierr); 2566f4f5c14SMatthew G. Knepley ierr = EG_getRange(edge, range, periodic); CHKERRQ(ierr); 2576f4f5c14SMatthew G. Knepley params[0] = 0.5*(range[0] + range[1]); 2586f4f5c14SMatthew G. Knepley ierr = EG_evaluate(edge, params, result);CHKERRQ(ierr); 2596f4f5c14SMatthew G. Knepley coords[v*cdim+0] = result[0]; 2606f4f5c14SMatthew G. Knepley coords[v*cdim+1] = result[1]; 2616f4f5c14SMatthew G. Knepley coords[v*cdim+2] = result[2]; 2626f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " Edge ID = %d (%D) \n", eid, v); 2636f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " (x,y,z) = (%lf, %lf, %lf)\n", coords[v*cdim+0], coords[v*cdim+1],coords[v*cdim+2]); 2646f4f5c14SMatthew G. Knepley params[0] = range[0]; 2656f4f5c14SMatthew G. Knepley ierr = EG_evaluate(edge, params, result);CHKERRQ(ierr); 2666f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " between (%lf, %lf, %lf)", result[0], result[1], result[2]); 2676f4f5c14SMatthew G. Knepley params[0] = range[1]; 2686f4f5c14SMatthew G. Knepley ierr = EG_evaluate(edge, params, result);CHKERRQ(ierr); 2696f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " and (%lf, %lf, %lf)\n\n", result[0], result[1], result[2]); 2706f4f5c14SMatthew G. Knepley } 2716f4f5c14SMatthew G. Knepley } 2726f4f5c14SMatthew G. Knepley if (Ner == 4) { 2736f4f5c14SMatthew G. Knepley PetscInt v = fOff + numQuads++; 2746f4f5c14SMatthew G. Knepley ego *fobjs; 2756f4f5c14SMatthew G. Knepley double range[4], params[3] = {0., 0., 0.}, result[18]; 2766f4f5c14SMatthew G. Knepley int Nf, face, periodic[2]; 2776f4f5c14SMatthew G. Knepley 2786f4f5c14SMatthew G. Knepley ierr = EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs);CHKERRQ(ierr); 2796f4f5c14SMatthew G. Knepley if (Nf != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Loop %d has %d faces, instead of 1", lid-1, Nf); 2806f4f5c14SMatthew G. Knepley face = EG_indexBodyTopo(body, fobjs[0]);CHKERRQ(ierr); 2816f4f5c14SMatthew G. Knepley ierr = EG_getRange(fobjs[0], range, periodic); CHKERRQ(ierr); 2826f4f5c14SMatthew G. Knepley params[0] = 0.5*(range[0] + range[1]); 2836f4f5c14SMatthew G. Knepley params[1] = 0.5*(range[2] + range[3]); 2846f4f5c14SMatthew G. Knepley ierr = EG_evaluate(fobjs[0], params, result);CHKERRQ(ierr); 2856f4f5c14SMatthew G. Knepley coords[v*cdim+0] = result[0]; 2866f4f5c14SMatthew G. Knepley coords[v*cdim+1] = result[1]; 2876f4f5c14SMatthew G. Knepley coords[v*cdim+2] = result[2]; 2886f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " Face ID = %d (%D) \n", face-1, v); 2896f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " (x,y,z) = (%lf, %lf, %lf) \n \n", coords[v*cdim+0], coords[v*cdim+1],coords[v*cdim+2]); 2906f4f5c14SMatthew G. Knepley } 2916f4f5c14SMatthew G. Knepley } 2926f4f5c14SMatthew G. Knepley } 2936f4f5c14SMatthew G. Knepley if (numEdges + numQuads != newVertices) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of new vertices %D != %D previous count", numEdges + numQuads, newVertices); 2946f4f5c14SMatthew G. Knepley CHKMEMQ; 295c4762a1bSJed Brown 296c4762a1bSJed Brown /* Get cell vertices by traversing loops */ 2976f4f5c14SMatthew G. Knepley numQuads = 0; 2986f4f5c14SMatthew G. Knepley cOff = 0; 299c4762a1bSJed Brown for (b = 0; b < nbodies; ++b) { 300c4762a1bSJed Brown ego body = bodies[b]; 301c4762a1bSJed Brown int id, Nl, l; 302c4762a1bSJed Brown 303c4762a1bSJed Brown ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr); 304c4762a1bSJed Brown for (l = 0; l < Nl; ++l) { 305c4762a1bSJed Brown ego loop = lobjs[l]; 3066f4f5c14SMatthew G. Knepley int lid, Ner = 0, Ne, e, nc = 0, c, Nt, t; 307c4762a1bSJed Brown 308c4762a1bSJed Brown lid = EG_indexBodyTopo(body, loop);CHKERRQ(ierr); 309c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, " LOOP ID: %d \n", lid);CHKERRQ(ierr); 310c4762a1bSJed Brown ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr); 311c4762a1bSJed Brown 312c4762a1bSJed Brown for (e = 0; e < Ne; ++e) { 313c4762a1bSJed Brown ego edge = objs[e]; 3146f4f5c14SMatthew G. Knepley int points[3]; 3156f4f5c14SMatthew G. Knepley int eid, Nv, v, tmp; 316c4762a1bSJed Brown 3176f4f5c14SMatthew G. Knepley eid = EG_indexBodyTopo(body, edge); 3186f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " EDGE ID: %d \n", eid);CHKERRQ(ierr); 319c4762a1bSJed Brown ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses); 3206f4f5c14SMatthew G. Knepley if (mtype == DEGENERATE) {ierr = PetscPrintf(PETSC_COMM_SELF, " EGDE %d is DEGENERATE \n", eid);CHKERRQ(ierr); continue;} 3216f4f5c14SMatthew G. Knepley else {++Ner;} 3226f4f5c14SMatthew G. Knepley if (Nv != 2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Edge %d has %d vertices != 2", eid, Nv); 323c4762a1bSJed Brown 324c4762a1bSJed Brown for (v = 0; v < Nv; ++v) { 325c4762a1bSJed Brown ego vertex = nobjs[v]; 326c4762a1bSJed Brown 327c4762a1bSJed Brown id = EG_indexBodyTopo(body, vertex); 3286f4f5c14SMatthew G. Knepley points[v*2] = id-1; 3296f4f5c14SMatthew G. Knepley } 3306f4f5c14SMatthew G. Knepley { 3316f4f5c14SMatthew G. Knepley PetscInt edgeNum; 3326f4f5c14SMatthew G. Knepley 3336f4f5c14SMatthew G. Knepley ierr = PetscHMapIGet(edgeMap, eid-1, &edgeNum);CHKERRQ(ierr); 3346f4f5c14SMatthew G. Knepley points[1] = numVertices - newVertices + edgeNum; 3356f4f5c14SMatthew G. Knepley } 3366f4f5c14SMatthew G. Knepley /* EGADS loops are not oriented, but seem to be in order, so we must piece them together */ 3376f4f5c14SMatthew G. Knepley if (!nc) { 3386f4f5c14SMatthew G. Knepley for (v = 0; v < Nv+1; ++v) cone[nc++] = points[v]; 3396f4f5c14SMatthew G. Knepley } else { 3406f4f5c14SMatthew G. Knepley if (cone[nc-1] == points[0]) {cone[nc++] = points[1]; if (cone[0] != points[2]) cone[nc++] = points[2];} 3416f4f5c14SMatthew G. Knepley else if (cone[nc-1] == points[2]) {cone[nc++] = points[1]; if (cone[0] != points[0]) cone[nc++] = points[0];} 3426f4f5c14SMatthew G. 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];} 3436f4f5c14SMatthew G. 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];} 3446f4f5c14SMatthew G. Knepley else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %d does not match its predecessor", eid); 345c4762a1bSJed Brown } 346c4762a1bSJed Brown } 3476f4f5c14SMatthew G. Knepley if (nc != 2*Ner) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %D != %D", nc, 2*Ner); 3486f4f5c14SMatthew G. Knepley if (Ner == 4) {cone[nc++] = numVertices - newVertices + numEdges + numQuads++;} 3496f4f5c14SMatthew G. Knepley if (nc > maxCorners) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Number of corners %D > %D max", nc, maxCorners); 3506f4f5c14SMatthew G. Knepley /* Triangulate the loop */ 3516f4f5c14SMatthew G. Knepley switch (Ner) { 3526f4f5c14SMatthew G. Knepley case 2: /* Bi-Segment -> 2 triangles */ 3536f4f5c14SMatthew G. Knepley Nt = 2; 3546f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[0]; 3556f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[1]; 3566f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[2]; 3576f4f5c14SMatthew G. Knepley ++cOff; 3586f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[0]; 3596f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[2]; 3606f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[3]; 3616f4f5c14SMatthew G. Knepley ++cOff; 3626f4f5c14SMatthew G. Knepley break; 3636f4f5c14SMatthew G. Knepley case 3: /* Triangle -> 4 triangles */ 3646f4f5c14SMatthew G. Knepley Nt = 4; 3656f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[0]; 3666f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[1]; 3676f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[5]; 3686f4f5c14SMatthew G. Knepley ++cOff; 3696f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[1]; 3706f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[2]; 3716f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[3]; 3726f4f5c14SMatthew G. Knepley ++cOff; 3736f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[5]; 3746f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[3]; 3756f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[4]; 3766f4f5c14SMatthew G. Knepley ++cOff; 3776f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[1]; 3786f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[3]; 3796f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[5]; 3806f4f5c14SMatthew G. Knepley ++cOff; 3816f4f5c14SMatthew G. Knepley break; 3826f4f5c14SMatthew G. Knepley case 4: /* Quad -> 8 triangles */ 3836f4f5c14SMatthew G. Knepley Nt = 8; 3846f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[0]; 3856f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[1]; 3866f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[7]; 3876f4f5c14SMatthew G. Knepley ++cOff; 3886f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[1]; 3896f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[2]; 3906f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[3]; 3916f4f5c14SMatthew G. Knepley ++cOff; 3926f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[3]; 3936f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[4]; 3946f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[5]; 3956f4f5c14SMatthew G. Knepley ++cOff; 3966f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[5]; 3976f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[6]; 3986f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[7]; 3996f4f5c14SMatthew G. Knepley ++cOff; 4006f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[8]; 4016f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[1]; 4026f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[3]; 4036f4f5c14SMatthew G. Knepley ++cOff; 4046f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[8]; 4056f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[3]; 4066f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[5]; 4076f4f5c14SMatthew G. Knepley ++cOff; 4086f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[8]; 4096f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[5]; 4106f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[7]; 4116f4f5c14SMatthew G. Knepley ++cOff; 4126f4f5c14SMatthew G. Knepley cells[cOff*numCorners+0] = cone[8]; 4136f4f5c14SMatthew G. Knepley cells[cOff*numCorners+1] = cone[7]; 4146f4f5c14SMatthew G. Knepley cells[cOff*numCorners+2] = cone[1]; 4156f4f5c14SMatthew G. Knepley ++cOff; 4166f4f5c14SMatthew G. Knepley break; 4176f4f5c14SMatthew G. Knepley default: SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d edges, which we do not support", lid, Ner); 4186f4f5c14SMatthew G. Knepley } 4196f4f5c14SMatthew G. Knepley for (t = 0; t < Nt; ++t) { 4206f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, " LOOP Corner NODEs Triangle %D (", t); 421c4762a1bSJed Brown for (c = 0; c < numCorners; ++c) { 422c4762a1bSJed Brown if (c > 0) {ierr = PetscPrintf(PETSC_COMM_SELF, ", ");} 4236f4f5c14SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_SELF, "%D", cells[(cOff-Nt+t)*numCorners+c]); 424c4762a1bSJed Brown } 425c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF, ")\n"); 426c4762a1bSJed Brown } 427c4762a1bSJed Brown } 4286f4f5c14SMatthew G. Knepley EG_free(lobjs); 429c4762a1bSJed Brown } 4306f4f5c14SMatthew G. Knepley } 4316f4f5c14SMatthew G. Knepley if (cOff != numCells) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Count of total cells %D != %D previous count", cOff, numCells); 432c4762a1bSJed Brown ierr = DMPlexCreateFromCellList(PETSC_COMM_WORLD, dim, numCells, numVertices, numCorners, PETSC_TRUE, cells, cdim, coords, &dm);CHKERRQ(ierr); 4336f4f5c14SMatthew G. Knepley ierr = PetscFree3(coords, cells, cone);CHKERRQ(ierr); 4346f4f5c14SMatthew G. Knepley /* Embed EGADS model in DM */ 435c4762a1bSJed Brown { 436c4762a1bSJed Brown PetscContainer modelObj; 437c4762a1bSJed Brown 438c4762a1bSJed Brown ierr = PetscContainerCreate(PETSC_COMM_SELF, &modelObj);CHKERRQ(ierr); 439c4762a1bSJed Brown ierr = PetscContainerSetPointer(modelObj, model);CHKERRQ(ierr); 440c4762a1bSJed Brown ierr = PetscObjectCompose((PetscObject) dm, "EGADS Model", (PetscObject) modelObj);CHKERRQ(ierr); 441c4762a1bSJed Brown ierr = PetscContainerDestroy(&modelObj);CHKERRQ(ierr); 442c4762a1bSJed Brown } 4436f4f5c14SMatthew G. Knepley /* Label points */ 444c4762a1bSJed Brown ierr = DMCreateLabel(dm, "EGADS Body ID");CHKERRQ(ierr); 445c4762a1bSJed Brown ierr = DMGetLabel(dm, "EGADS Body ID", &bodyLabel);CHKERRQ(ierr); 446c4762a1bSJed Brown ierr = DMCreateLabel(dm, "EGADS Face ID");CHKERRQ(ierr); 447c4762a1bSJed Brown ierr = DMGetLabel(dm, "EGADS Face ID", &faceLabel);CHKERRQ(ierr); 448c4762a1bSJed Brown ierr = DMCreateLabel(dm, "EGADS Edge ID");CHKERRQ(ierr); 449c4762a1bSJed Brown ierr = DMGetLabel(dm, "EGADS Edge ID", &edgeLabel);CHKERRQ(ierr); 4506f4f5c14SMatthew G. Knepley cOff = 0; 451c4762a1bSJed Brown for (b = 0; b < nbodies; ++b) { 452c4762a1bSJed Brown ego body = bodies[b]; 453c4762a1bSJed Brown int id, Nl, l; 454c4762a1bSJed Brown 455c4762a1bSJed Brown ierr = EG_getBodyTopos(body, NULL, LOOP, &Nl, &lobjs);CHKERRQ(ierr); 456c4762a1bSJed Brown for (l = 0; l < Nl; ++l) { 457c4762a1bSJed Brown ego loop = lobjs[l]; 4586f4f5c14SMatthew G. Knepley ego *fobjs; 4596f4f5c14SMatthew G. Knepley int lid, Nf, fid, Ner = 0, Ne, e, Nt = 0, t; 460c4762a1bSJed Brown 461c4762a1bSJed Brown lid = EG_indexBodyTopo(body, loop);CHKERRQ(ierr); 462c4762a1bSJed Brown ierr = EG_getBodyTopos(body, loop, FACE, &Nf, &fobjs);CHKERRQ(ierr); 4636f4f5c14SMatthew G. Knepley if (Nf > 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Loop %d has %d > 1 faces, which is not supported", lid, Nf);CHKERRQ(ierr); 464c4762a1bSJed Brown fid = EG_indexBodyTopo(body, fobjs[0]);CHKERRQ(ierr); 4656f4f5c14SMatthew G. Knepley EG_free(fobjs); 466c4762a1bSJed Brown ierr = EG_getTopology(loop, &geom, &oclass, &mtype, NULL, &Ne, &objs, &senses);CHKERRQ(ierr); 467c4762a1bSJed Brown for (e = 0; e < Ne; ++e) { 468c4762a1bSJed Brown ego edge = objs[e]; 469c4762a1bSJed Brown int eid, Nv, v; 4706f4f5c14SMatthew G. Knepley PetscInt points[3], support[2], numEdges, edgeNum; 471c4762a1bSJed Brown const PetscInt *edges; 472c4762a1bSJed Brown 473c4762a1bSJed Brown eid = EG_indexBodyTopo(body, edge); 474c4762a1bSJed Brown ierr = EG_getTopology(edge, &geom, &oclass, &mtype, NULL, &Nv, &nobjs, &senses); 4756f4f5c14SMatthew G. Knepley if (mtype == DEGENERATE) {continue;} 4766f4f5c14SMatthew G. Knepley else {++Ner;} 477c4762a1bSJed Brown for (v = 0; v < Nv; ++v) { 478c4762a1bSJed Brown ego vertex = nobjs[v]; 479c4762a1bSJed Brown 480c4762a1bSJed Brown id = EG_indexBodyTopo(body, vertex); 481c4762a1bSJed Brown ierr = DMLabelSetValue(edgeLabel, numCells + id-1, eid);CHKERRQ(ierr); 4826f4f5c14SMatthew G. Knepley points[v*2] = numCells + id-1; 483c4762a1bSJed Brown } 4846f4f5c14SMatthew G. Knepley ierr = PetscHMapIGet(edgeMap, eid-1, &edgeNum);CHKERRQ(ierr); 4856f4f5c14SMatthew G. Knepley points[1] = numCells + numVertices - newVertices + edgeNum; 4866f4f5c14SMatthew G. Knepley 4876f4f5c14SMatthew G. Knepley ierr = DMLabelSetValue(edgeLabel, points[1], eid);CHKERRQ(ierr); 4886f4f5c14SMatthew G. Knepley support[0] = points[0]; 4896f4f5c14SMatthew G. Knepley support[1] = points[1]; 490c4762a1bSJed Brown ierr = DMPlexGetJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr); 4916f4f5c14SMatthew G. Knepley if (numEdges != 1) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vertices (%D, %D) should only bound 1 edge, not %D", support[0], support[1], numEdges); 4926f4f5c14SMatthew G. Knepley ierr = DMLabelSetValue(edgeLabel, edges[0], eid);CHKERRQ(ierr); 4936f4f5c14SMatthew G. Knepley ierr = DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr); 4946f4f5c14SMatthew G. Knepley support[0] = points[1]; 4956f4f5c14SMatthew G. Knepley support[1] = points[2]; 4966f4f5c14SMatthew G. Knepley ierr = DMPlexGetJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr); 4976f4f5c14SMatthew G. Knepley if (numEdges != 1) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vertices (%D, %D) should only bound 1 edge, not %D", support[0], support[1], numEdges); 498c4762a1bSJed Brown ierr = DMLabelSetValue(edgeLabel, edges[0], eid);CHKERRQ(ierr); 499c4762a1bSJed Brown ierr = DMPlexRestoreJoin(dm, 2, support, &numEdges, &edges);CHKERRQ(ierr); 500c4762a1bSJed Brown } 5016f4f5c14SMatthew G. Knepley switch (Ner) { 5026f4f5c14SMatthew G. Knepley case 2: Nt = 2;break; 5036f4f5c14SMatthew G. Knepley case 3: Nt = 4;break; 5046f4f5c14SMatthew G. Knepley case 4: Nt = 8;break; 5056f4f5c14SMatthew G. Knepley default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Loop with %d edges is unsupported", Ner); 506c4762a1bSJed Brown } 5076f4f5c14SMatthew G. Knepley for (t = 0; t < Nt; ++t) { 5086f4f5c14SMatthew G. Knepley ierr = DMLabelSetValue(bodyLabel, cOff+t, b);CHKERRQ(ierr); 5096f4f5c14SMatthew G. Knepley ierr = DMLabelSetValue(faceLabel, cOff+t, fid);CHKERRQ(ierr); 510c4762a1bSJed Brown } 5116f4f5c14SMatthew G. Knepley cOff += Nt; 512c4762a1bSJed Brown } 5136f4f5c14SMatthew G. Knepley EG_free(lobjs); 5146f4f5c14SMatthew G. Knepley } 5156f4f5c14SMatthew G. Knepley ierr = PetscHMapIDestroy(&edgeMap);CHKERRQ(ierr); 516c4762a1bSJed Brown ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 517c4762a1bSJed Brown for (c = cStart; c < cEnd; ++c) { 518c4762a1bSJed Brown PetscInt *closure = NULL; 519c4762a1bSJed Brown PetscInt clSize, cl, bval, fval; 520c4762a1bSJed Brown 521c4762a1bSJed Brown ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr); 522c4762a1bSJed Brown ierr = DMLabelGetValue(bodyLabel, c, &bval);CHKERRQ(ierr); 523c4762a1bSJed Brown ierr = DMLabelGetValue(faceLabel, c, &fval);CHKERRQ(ierr); 524c4762a1bSJed Brown for (cl = 0; cl < clSize*2; cl += 2) { 525c4762a1bSJed Brown ierr = DMLabelSetValue(bodyLabel, closure[cl], bval);CHKERRQ(ierr); 526c4762a1bSJed Brown ierr = DMLabelSetValue(faceLabel, closure[cl], fval);CHKERRQ(ierr); 527c4762a1bSJed Brown } 528c4762a1bSJed Brown ierr = DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &clSize, &closure);CHKERRQ(ierr); 529c4762a1bSJed Brown } 530c4762a1bSJed Brown ierr = DMLabelView(bodyLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 531c4762a1bSJed Brown ierr = DMLabelView(faceLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 532c4762a1bSJed Brown ierr = DMLabelView(edgeLabel, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 533c4762a1bSJed Brown ierr = DMSetFromOptions(dm);CHKERRQ(ierr); 534c4762a1bSJed Brown 535c4762a1bSJed Brown ierr = DMViewFromOptions(dm, NULL, "-dm_view");CHKERRQ(ierr); 536c4762a1bSJed Brown ierr = DMDestroy(&dm);CHKERRQ(ierr); 537c4762a1bSJed Brown 538c4762a1bSJed Brown /* Close EGADSlite file */ 539c4762a1bSJed Brown ierr = EG_close(context);CHKERRQ(ierr); 540c4762a1bSJed Brown ierr = PetscFinalize(); 541c4762a1bSJed Brown return ierr; 542c4762a1bSJed Brown } 543c4762a1bSJed Brown 544c4762a1bSJed Brown /*TEST 545c4762a1bSJed Brown 546c4762a1bSJed Brown build: 547c4762a1bSJed Brown requires: egads 548c4762a1bSJed Brown 549c4762a1bSJed Brown test: 550c4762a1bSJed Brown suffix: sphere_0 5516f4f5c14SMatthew G. Knepley args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/unit_sphere.egadslite -dm_refine 1 -dm_view -dm_plex_check_all 5526f4f5c14SMatthew G. Knepley 5536f4f5c14SMatthew G. Knepley test: 5546f4f5c14SMatthew G. Knepley suffix: nozzle_0 5556f4f5c14SMatthew G. Knepley args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/nozzle.egadslite -dm_refine 1 -dm_view -dm_plex_check_all 556c4762a1bSJed Brown 557c4762a1bSJed Brown TEST*/ 558