1*2e1d0745SJose E. Roman #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*2e1d0745SJose E. Roman #include <petsc/private/isimpl.h> 3*2e1d0745SJose E. Roman #include <petsc/private/vecimpl.h> 4*2e1d0745SJose E. Roman #include <petsclayouthdf5.h> 5*2e1d0745SJose E. Roman 6*2e1d0745SJose E. Roman static PetscErrorCode SplitPath_Private(char path[], char name[]) 7*2e1d0745SJose E. Roman { 8*2e1d0745SJose E. Roman char *tmp; 9*2e1d0745SJose E. Roman size_t len; 10*2e1d0745SJose E. Roman 11*2e1d0745SJose E. Roman PetscFunctionBegin; 12*2e1d0745SJose E. Roman PetscCall(PetscStrrchr(path, '/', &tmp)); 13*2e1d0745SJose E. Roman PetscCall(PetscStrlen(tmp, &len)); 14*2e1d0745SJose E. Roman PetscCall(PetscStrncpy(name, tmp, len + 1)); /* assuming adequate buffer */ 15*2e1d0745SJose E. Roman if (tmp != path) { 16*2e1d0745SJose E. Roman /* '/' found, name is substring of path after last occurrence of '/'. */ 17*2e1d0745SJose E. Roman /* Trim the '/name' part from path just by inserting null character. */ 18*2e1d0745SJose E. Roman tmp--; 19*2e1d0745SJose E. Roman *tmp = '\0'; 20*2e1d0745SJose E. Roman } else { 21*2e1d0745SJose E. Roman /* '/' not found, name = path, path = "/". */ 22*2e1d0745SJose E. Roman PetscCall(PetscStrncpy(path, "/", 2)); /* assuming adequate buffer */ 23*2e1d0745SJose E. Roman } 24*2e1d0745SJose E. Roman PetscFunctionReturn(PETSC_SUCCESS); 25*2e1d0745SJose E. Roman } 26*2e1d0745SJose E. Roman 27*2e1d0745SJose E. Roman /* 28*2e1d0745SJose E. Roman - invert (involute) cells of some types according to XDMF/VTK numbering of vertices in a cells 29*2e1d0745SJose E. Roman - cell type is identified using the number of vertices 30*2e1d0745SJose E. Roman */ 31*2e1d0745SJose E. Roman static PetscErrorCode DMPlexInvertCells_XDMF_Private(DM dm) 32*2e1d0745SJose E. Roman { 33*2e1d0745SJose E. Roman PetscInt dim, *cones, cHeight, cStart, cEnd, p; 34*2e1d0745SJose E. Roman PetscSection cs; 35*2e1d0745SJose E. Roman 36*2e1d0745SJose E. Roman PetscFunctionBegin; 37*2e1d0745SJose E. Roman PetscCall(DMGetDimension(dm, &dim)); 38*2e1d0745SJose E. Roman if (dim != 3) PetscFunctionReturn(PETSC_SUCCESS); 39*2e1d0745SJose E. Roman PetscCall(DMPlexGetCones(dm, &cones)); 40*2e1d0745SJose E. Roman PetscCall(DMPlexGetConeSection(dm, &cs)); 41*2e1d0745SJose E. Roman PetscCall(DMPlexGetVTKCellHeight(dm, &cHeight)); 42*2e1d0745SJose E. Roman PetscCall(DMPlexGetHeightStratum(dm, cHeight, &cStart, &cEnd)); 43*2e1d0745SJose E. Roman for (p = cStart; p < cEnd; p++) { 44*2e1d0745SJose E. Roman PetscInt numCorners, o; 45*2e1d0745SJose E. Roman 46*2e1d0745SJose E. Roman PetscCall(PetscSectionGetDof(cs, p, &numCorners)); 47*2e1d0745SJose E. Roman PetscCall(PetscSectionGetOffset(cs, p, &o)); 48*2e1d0745SJose E. Roman switch (numCorners) { 49*2e1d0745SJose E. Roman case 4: 50*2e1d0745SJose E. Roman PetscCall(DMPlexInvertCell(DM_POLYTOPE_TETRAHEDRON, &cones[o])); 51*2e1d0745SJose E. Roman break; 52*2e1d0745SJose E. Roman case 6: 53*2e1d0745SJose E. Roman PetscCall(DMPlexInvertCell(DM_POLYTOPE_TRI_PRISM, &cones[o])); 54*2e1d0745SJose E. Roman break; 55*2e1d0745SJose E. Roman case 8: 56*2e1d0745SJose E. Roman PetscCall(DMPlexInvertCell(DM_POLYTOPE_HEXAHEDRON, &cones[o])); 57*2e1d0745SJose E. Roman break; 58*2e1d0745SJose E. Roman } 59*2e1d0745SJose E. Roman } 60*2e1d0745SJose E. Roman PetscFunctionReturn(PETSC_SUCCESS); 61*2e1d0745SJose E. Roman } 62*2e1d0745SJose E. Roman 63*2e1d0745SJose E. Roman PetscErrorCode DMPlexLoad_HDF5_Xdmf_Internal(DM dm, PetscViewer viewer) 64*2e1d0745SJose E. Roman { 65*2e1d0745SJose E. Roman Vec coordinates; 66*2e1d0745SJose E. Roman IS cells; 67*2e1d0745SJose E. Roman PetscInt spatialDim, topoDim = -1, numCells, numVertices, NVertices, numCorners; 68*2e1d0745SJose E. Roman PetscMPIInt rank; 69*2e1d0745SJose E. Roman MPI_Comm comm; 70*2e1d0745SJose E. Roman char topo_path[PETSC_MAX_PATH_LEN] = "/viz/topology/cells", topo_name[PETSC_MAX_PATH_LEN]; 71*2e1d0745SJose E. Roman char geom_path[PETSC_MAX_PATH_LEN] = "/geometry/vertices", geom_name[PETSC_MAX_PATH_LEN]; 72*2e1d0745SJose E. Roman PetscBool seq = PETSC_FALSE; 73*2e1d0745SJose E. Roman 74*2e1d0745SJose E. Roman PetscFunctionBegin; 75*2e1d0745SJose E. Roman PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 76*2e1d0745SJose E. Roman PetscCallMPI(MPI_Comm_rank(comm, &rank)); 77*2e1d0745SJose E. Roman 78*2e1d0745SJose E. Roman PetscOptionsBegin(PetscObjectComm((PetscObject)dm), ((PetscObject)dm)->prefix, "DMPlex HDF5/XDMF Loader Options", "PetscViewer"); 79*2e1d0745SJose E. Roman PetscCall(PetscOptionsString("-dm_plex_hdf5_topology_path", "HDF5 path of topology dataset", NULL, topo_path, topo_path, sizeof(topo_path), NULL)); 80*2e1d0745SJose E. Roman PetscCall(PetscOptionsString("-dm_plex_hdf5_geometry_path", "HDF5 path to geometry dataset", NULL, geom_path, geom_path, sizeof(geom_path), NULL)); 81*2e1d0745SJose E. Roman PetscCall(PetscOptionsBool("-dm_plex_hdf5_force_sequential", "force sequential loading", NULL, seq, &seq, NULL)); 82*2e1d0745SJose E. Roman PetscOptionsEnd(); 83*2e1d0745SJose E. Roman 84*2e1d0745SJose E. Roman PetscCall(SplitPath_Private(topo_path, topo_name)); 85*2e1d0745SJose E. Roman PetscCall(SplitPath_Private(geom_path, geom_name)); 86*2e1d0745SJose E. Roman PetscCall(PetscInfo(dm, "Topology group %s, name %s\n", topo_path, topo_name)); 87*2e1d0745SJose E. Roman PetscCall(PetscInfo(dm, "Geometry group %s, name %s\n", geom_path, geom_name)); 88*2e1d0745SJose E. Roman 89*2e1d0745SJose E. Roman /* Read topology */ 90*2e1d0745SJose E. Roman PetscCall(PetscViewerHDF5PushGroup(viewer, topo_path)); 91*2e1d0745SJose E. Roman PetscCall(ISCreate(comm, &cells)); 92*2e1d0745SJose E. Roman PetscCall(PetscObjectSetName((PetscObject)cells, topo_name)); 93*2e1d0745SJose E. Roman if (seq) { 94*2e1d0745SJose E. Roman PetscCall(PetscViewerHDF5ReadSizes(viewer, topo_name, NULL, &numCells)); 95*2e1d0745SJose E. Roman PetscCall(PetscLayoutSetSize(cells->map, numCells)); 96*2e1d0745SJose E. Roman numCells = rank == 0 ? numCells : 0; 97*2e1d0745SJose E. Roman PetscCall(PetscLayoutSetLocalSize(cells->map, numCells)); 98*2e1d0745SJose E. Roman } 99*2e1d0745SJose E. Roman PetscCall(ISLoad(cells, viewer)); 100*2e1d0745SJose E. Roman PetscCall(ISGetLocalSize(cells, &numCells)); 101*2e1d0745SJose E. Roman PetscCall(ISGetBlockSize(cells, &numCorners)); 102*2e1d0745SJose E. Roman PetscCall(PetscViewerHDF5ReadAttribute(viewer, topo_name, "cell_dim", PETSC_INT, &topoDim, &topoDim)); 103*2e1d0745SJose E. Roman PetscCall(PetscViewerHDF5PopGroup(viewer)); 104*2e1d0745SJose E. Roman numCells /= numCorners; 105*2e1d0745SJose E. Roman 106*2e1d0745SJose E. Roman /* Read geometry */ 107*2e1d0745SJose E. Roman PetscCall(PetscViewerHDF5PushGroup(viewer, geom_path)); 108*2e1d0745SJose E. Roman PetscCall(VecCreate(comm, &coordinates)); 109*2e1d0745SJose E. Roman PetscCall(PetscObjectSetName((PetscObject)coordinates, geom_name)); 110*2e1d0745SJose E. Roman if (seq) { 111*2e1d0745SJose E. Roman PetscCall(PetscViewerHDF5ReadSizes(viewer, geom_name, NULL, &numVertices)); 112*2e1d0745SJose E. Roman PetscCall(PetscLayoutSetSize(coordinates->map, numVertices)); 113*2e1d0745SJose E. Roman numVertices = rank == 0 ? numVertices : 0; 114*2e1d0745SJose E. Roman PetscCall(PetscLayoutSetLocalSize(coordinates->map, numVertices)); 115*2e1d0745SJose E. Roman } 116*2e1d0745SJose E. Roman PetscCall(VecLoad(coordinates, viewer)); 117*2e1d0745SJose E. Roman PetscCall(VecGetLocalSize(coordinates, &numVertices)); 118*2e1d0745SJose E. Roman PetscCall(VecGetSize(coordinates, &NVertices)); 119*2e1d0745SJose E. Roman PetscCall(VecGetBlockSize(coordinates, &spatialDim)); 120*2e1d0745SJose E. Roman PetscCall(PetscViewerHDF5PopGroup(viewer)); 121*2e1d0745SJose E. Roman numVertices /= spatialDim; 122*2e1d0745SJose E. Roman NVertices /= spatialDim; 123*2e1d0745SJose E. Roman 124*2e1d0745SJose E. Roman PetscCall(PetscInfo(NULL, "Loaded mesh dimensions: numCells %" PetscInt_FMT " numCorners %" PetscInt_FMT " numVertices %" PetscInt_FMT " spatialDim %" PetscInt_FMT "\n", numCells, numCorners, numVertices, spatialDim)); 125*2e1d0745SJose E. Roman { 126*2e1d0745SJose E. Roman const PetscScalar *coordinates_arr; 127*2e1d0745SJose E. Roman PetscReal *coordinates_arr_real; 128*2e1d0745SJose E. Roman const PetscInt *cells_arr; 129*2e1d0745SJose E. Roman PetscSF sfVert = NULL; 130*2e1d0745SJose E. Roman PetscInt i; 131*2e1d0745SJose E. Roman 132*2e1d0745SJose E. Roman PetscCall(VecGetArrayRead(coordinates, &coordinates_arr)); 133*2e1d0745SJose E. Roman PetscCall(ISGetIndices(cells, &cells_arr)); 134*2e1d0745SJose E. Roman 135*2e1d0745SJose E. Roman if (PetscDefined(USE_COMPLEX)) { 136*2e1d0745SJose E. Roman /* convert to real numbers if PetscScalar is complex */ 137*2e1d0745SJose E. Roman /*TODO More systematic would be to change all the function arguments to PetscScalar */ 138*2e1d0745SJose E. Roman PetscCall(PetscMalloc1(numVertices * spatialDim, &coordinates_arr_real)); 139*2e1d0745SJose E. Roman for (i = 0; i < numVertices * spatialDim; ++i) { 140*2e1d0745SJose E. Roman coordinates_arr_real[i] = PetscRealPart(coordinates_arr[i]); 141*2e1d0745SJose E. Roman PetscAssert(PetscImaginaryPart(coordinates_arr[i]) == 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Vector of coordinates contains complex numbers but only real vectors are currently supported."); 142*2e1d0745SJose E. Roman } 143*2e1d0745SJose E. Roman } else coordinates_arr_real = (PetscReal *)coordinates_arr; 144*2e1d0745SJose E. Roman 145*2e1d0745SJose E. Roman PetscCall(DMSetDimension(dm, topoDim < 0 ? spatialDim : topoDim)); 146*2e1d0745SJose E. Roman PetscCall(DMPlexBuildFromCellListParallel(dm, numCells, numVertices, NVertices, numCorners, cells_arr, &sfVert, NULL)); 147*2e1d0745SJose E. Roman PetscCall(DMPlexInvertCells_XDMF_Private(dm)); 148*2e1d0745SJose E. Roman PetscCall(DMPlexBuildCoordinatesFromCellListParallel(dm, spatialDim, sfVert, coordinates_arr_real)); 149*2e1d0745SJose E. Roman PetscCall(VecRestoreArrayRead(coordinates, &coordinates_arr)); 150*2e1d0745SJose E. Roman PetscCall(ISRestoreIndices(cells, &cells_arr)); 151*2e1d0745SJose E. Roman PetscCall(PetscSFDestroy(&sfVert)); 152*2e1d0745SJose E. Roman if (PetscDefined(USE_COMPLEX)) PetscCall(PetscFree(coordinates_arr_real)); 153*2e1d0745SJose E. Roman } 154*2e1d0745SJose E. Roman PetscCall(ISDestroy(&cells)); 155*2e1d0745SJose E. Roman PetscCall(VecDestroy(&coordinates)); 156*2e1d0745SJose E. Roman 157*2e1d0745SJose E. Roman /* scale coordinates - unlike in DMPlexLoad_HDF5_Internal, this can only be done after DM is populated */ 158*2e1d0745SJose E. Roman { 159*2e1d0745SJose E. Roman PetscReal lengthScale; 160*2e1d0745SJose E. Roman 161*2e1d0745SJose E. Roman PetscCall(DMPlexGetScale(dm, PETSC_UNIT_LENGTH, &lengthScale)); 162*2e1d0745SJose E. Roman PetscCall(DMGetCoordinates(dm, &coordinates)); 163*2e1d0745SJose E. Roman PetscCall(VecScale(coordinates, 1.0 / lengthScale)); 164*2e1d0745SJose E. Roman } 165*2e1d0745SJose E. Roman 166*2e1d0745SJose E. Roman /* Read Labels */ 167*2e1d0745SJose E. Roman /* TODO: this probably does not work as elements get permuted */ 168*2e1d0745SJose E. Roman /* PetscCall(DMPlexLabelsLoad_HDF5_Internal(dm, viewer)); */ 169*2e1d0745SJose E. Roman PetscFunctionReturn(PETSC_SUCCESS); 170*2e1d0745SJose E. Roman } 171