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